xref: /petsc/src/tao/constrained/impls/almm/almm.c (revision 8fb5bd83c3955fefcf33a54e3bb66920a9fa884b)
1 #include <../src/tao/constrained/impls/almm/almm.h> /*I "petsctao.h" I*/ /*I "petscvec.h" I*/
2 #include <petsctao.h>
3 #include <petsc/private/petscimpl.h>
4 #include <petsc/private/vecimpl.h>
5 
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);
12 
13 static PetscErrorCode TaoSolve_ALMM(Tao tao)
14 {
15   TAO_ALMM           *auglag = (TAO_ALMM*)tao->data;
16   TaoConvergedReason reason;
17   PetscReal          updated;
18 
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(VecZeroEntries(auglag->Yi));
26     }
27     if (tao->eq_constrained) PetscCall(VecZeroEntries(auglag->Ye));
28   }
29 
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   PetscCall((*tao->ops->convergencetest)(tao, 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) {
46         PetscCall(VecDot(auglag->Ce, auglag->Ce, &auglag->cenorm));
47       }
48       auglag->cinorm = 0.0;
49       if (tao->ineq_constrained) {
50         PetscCall(VecCopy(auglag->Ci, auglag->Ciwork));
51         PetscCall(VecScale(auglag->Ciwork, -1.0));
52         PetscCall(VecPointwiseMax(auglag->Ciwork, auglag->Cizero, auglag->Ciwork));
53         PetscCall(VecDot(auglag->Ciwork, auglag->Ciwork, &auglag->cinorm));
54       }
55       /* determine initial penalty factor based on the balance of constraint violation and objective function value */
56       auglag->mu = PetscMax(1.e-6, PetscMin(10.0, 2.0*PetscAbsReal(auglag->fval)/(auglag->cenorm + auglag->cinorm)));
57       break;
58     default:
59       break;
60   }
61   auglag->gtol = auglag->gtol0;
62   PetscCall(PetscInfo(tao,"Initial penalty: %.2g\n",(double)auglag->mu));
63 
64   /* start aug-lag outer loop */
65   while (tao->reason == TAO_CONTINUE_ITERATING) {
66     ++tao->niter;
67     /* update subsolver tolerance */
68     PetscCall(PetscInfo(tao,"Subsolver tolerance: ||G|| <= %e\n",(double)auglag->gtol));
69     PetscCall(TaoSetTolerances(auglag->subsolver, auglag->gtol, 0.0, 0.0));
70     /* solve the bound-constrained or unconstrained subproblem */
71     PetscCall(TaoSolve(auglag->subsolver));
72     PetscCall(TaoGetConvergedReason(auglag->subsolver, &reason));
73     tao->ksp_its += auglag->subsolver->ksp_its;
74     if (reason != TAO_CONVERGED_GATOL) {
75       PetscCall(PetscInfo(tao,"Subsolver failed to converge, reason: %s\n",TaoConvergedReasons[reason]));
76     }
77     /* evaluate solution and test convergence */
78     PetscCall((*auglag->sub_obj)(tao));
79     PetscCall(TaoALMMComputeOptimalityNorms_Private(tao));
80     /* decide whether to update multipliers or not */
81     updated = 0.0;
82     if (auglag->cnorm <= auglag->ytol) {
83       PetscCall(PetscInfo(tao,"Multipliers updated: ||C|| <= %e\n",(double)auglag->ytol));
84       /* constraints are good, update multipliers and convergence tolerances */
85       if (tao->eq_constrained) {
86         PetscCall(VecAXPY(auglag->Ye, auglag->mu, auglag->Ce));
87         PetscCall(VecSet(auglag->Cework, auglag->ye_max));
88         PetscCall(VecPointwiseMin(auglag->Ye, auglag->Cework, auglag->Ye));
89         PetscCall(VecSet(auglag->Cework, auglag->ye_min));
90         PetscCall(VecPointwiseMax(auglag->Ye, auglag->Cework, auglag->Ye));
91       }
92       if (tao->ineq_constrained) {
93         PetscCall(VecAXPY(auglag->Yi, auglag->mu, auglag->Ci));
94         PetscCall(VecSet(auglag->Ciwork, auglag->yi_max));
95         PetscCall(VecPointwiseMin(auglag->Yi, auglag->Ciwork, auglag->Yi));
96         PetscCall(VecSet(auglag->Ciwork, auglag->yi_min));
97         PetscCall(VecPointwiseMax(auglag->Yi, auglag->Ciwork, auglag->Yi));
98       }
99       /* tolerances are updated only for non-PHR methods */
100       if (auglag->type != TAO_ALMM_PHR) {
101         auglag->ytol = PetscMax(tao->catol, auglag->ytol/PetscPowReal(auglag->mu, auglag->mu_pow_good));
102         auglag->gtol = PetscMax(tao->gatol, auglag->gtol/auglag->mu);
103       }
104       updated = 1.0;
105     } else {
106       /* constraints are bad, update penalty factor */
107       auglag->mu = PetscMin(auglag->mu_max, auglag->mu_fac*auglag->mu);
108       /* tolerances are reset only for non-PHR methods */
109       if (auglag->type != TAO_ALMM_PHR) {
110         auglag->ytol = PetscMax(tao->catol, 0.1/PetscPowReal(auglag->mu, auglag->mu_pow_bad));
111         auglag->gtol = PetscMax(tao->gatol, 1.0/auglag->mu);
112       }
113       PetscCall(PetscInfo(tao,"Penalty increased: mu = %.2g\n",(double)auglag->mu));
114     }
115     PetscCall(TaoLogConvergenceHistory(tao, auglag->fval, auglag->gnorm, auglag->cnorm, tao->ksp_its));
116     PetscCall(TaoMonitor(tao, tao->niter, auglag->fval, auglag->gnorm, auglag->cnorm, updated));
117     PetscCall((*tao->ops->convergencetest)(tao, tao->cnvP));
118   }
119 
120   PetscFunctionReturn(0);
121 }
122 
123 static PetscErrorCode TaoView_ALMM(Tao tao,PetscViewer viewer)
124 {
125   TAO_ALMM       *auglag = (TAO_ALMM*)tao->data;
126   PetscBool      isascii;
127 
128   PetscFunctionBegin;
129   PetscCall(PetscViewerASCIIPushTab(viewer));
130   PetscCall(TaoView(auglag->subsolver,viewer));
131   PetscCall(PetscViewerASCIIPopTab(viewer));
132   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
133   if (isascii) {
134     PetscCall(PetscViewerASCIIPushTab(viewer));
135     PetscCall(PetscViewerASCIIPrintf(viewer, "ALMM Formulation Type: %s\n", TaoALMMTypes[auglag->type]));
136     PetscCall(PetscViewerASCIIPopTab(viewer));
137   }
138   PetscFunctionReturn(0);
139 }
140 
141 static PetscErrorCode TaoSetUp_ALMM(Tao tao)
142 {
143   TAO_ALMM       *auglag = (TAO_ALMM*)tao->data;
144   VecType        vec_type;
145   Vec            SL, SU;
146   PetscBool      is_cg = PETSC_FALSE, is_lmvm = PETSC_FALSE;
147 
148   PetscFunctionBegin;
149   PetscCheck(!tao->ineq_doublesided,PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "TAOALMM does not support double-sided inequality constraint definition. Please restructure your inequality constrainst to fit the form c(x) >= 0.");
150   PetscCheck(tao->eq_constrained || tao->ineq_constrained,PetscObjectComm((PetscObject)tao), PETSC_ERR_ORDER, "Equality and/or inequality constraints must be defined before solver setup.");
151   PetscCall(TaoComputeVariableBounds(tao));
152   /* alias base vectors and create extras */
153   PetscCall(VecGetType(tao->solution, &vec_type));
154   auglag->Px = tao->solution;
155   if (!tao->gradient) { /* base gradient */
156     PetscCall(VecDuplicate(tao->solution, &tao->gradient));
157   }
158   auglag->LgradX = tao->gradient;
159   if (!auglag->Xwork) { /* opt var work vector */
160     PetscCall(VecDuplicate(tao->solution, &auglag->Xwork));
161   }
162   if (tao->eq_constrained) {
163     auglag->Ce = tao->constraints_equality;
164     auglag->Ae = tao->jacobian_equality;
165     if (!auglag->Ye) { /* equality multipliers */
166       PetscCall(VecDuplicate(auglag->Ce, &auglag->Ye));
167     }
168     if (!auglag->Cework) {
169       PetscCall(VecDuplicate(auglag->Ce, &auglag->Cework));
170     }
171   }
172   if (tao->ineq_constrained) {
173     auglag->Ci = tao->constraints_inequality;
174     auglag->Ai = tao->jacobian_inequality;
175     if (!auglag->Yi) { /* inequality multipliers */
176       PetscCall(VecDuplicate(auglag->Ci, &auglag->Yi));
177     }
178     if (!auglag->Ciwork) {
179       PetscCall(VecDuplicate(auglag->Ci, &auglag->Ciwork));
180     }
181     if (!auglag->Cizero) {
182       PetscCall(VecDuplicate(auglag->Ci, &auglag->Cizero));
183       PetscCall(VecZeroEntries(auglag->Cizero));
184     }
185     if (!auglag->Ps) { /* slack vars */
186       PetscCall(VecDuplicate(auglag->Ci, &auglag->Ps));
187     }
188     if (!auglag->LgradS) { /* slack component of Lagrangian gradient */
189       PetscCall(VecDuplicate(auglag->Ci, &auglag->LgradS));
190     }
191     /* create vector for combined primal space and the associated communication objects */
192     if (!auglag->P) {
193       PetscCall(PetscMalloc1(2, &auglag->Parr));
194       auglag->Parr[0] = auglag->Px; auglag->Parr[1] = auglag->Ps;
195       PetscCall(VecConcatenate(2, auglag->Parr, &auglag->P, &auglag->Pis));
196       PetscCall(PetscMalloc1(2, &auglag->Pscatter));
197       PetscCall(VecScatterCreate(auglag->P, auglag->Pis[0], auglag->Px, NULL, &auglag->Pscatter[0]));
198       PetscCall(VecScatterCreate(auglag->P, auglag->Pis[1], auglag->Ps, NULL, &auglag->Pscatter[1]));
199     }
200     if (tao->eq_constrained) {
201       /* create vector for combined dual space and the associated communication objects */
202       if (!auglag->Y) {
203         PetscCall(PetscMalloc1(2, &auglag->Yarr));
204         auglag->Yarr[0] = auglag->Ye; auglag->Yarr[1] = auglag->Yi;
205         PetscCall(VecConcatenate(2, auglag->Yarr, &auglag->Y, &auglag->Yis));
206         PetscCall(PetscMalloc1(2, &auglag->Yscatter));
207         PetscCall(VecScatterCreate(auglag->Y, auglag->Yis[0], auglag->Ye, NULL, &auglag->Yscatter[0]));
208         PetscCall(VecScatterCreate(auglag->Y, auglag->Yis[1], auglag->Yi, NULL, &auglag->Yscatter[1]));
209       }
210       if (!auglag->C) {
211         PetscCall(VecDuplicate(auglag->Y, &auglag->C));
212       }
213     } else {
214       if (!auglag->C) {
215         auglag->C = auglag->Ci;
216       }
217       if (!auglag->Y) {
218         auglag->Y = auglag->Yi;
219       }
220     }
221   } else {
222     if (!auglag->P) {
223       auglag->P = auglag->Px;
224     }
225     if (!auglag->G) {
226       auglag->G = auglag->LgradX;
227     }
228     if (!auglag->C) {
229       auglag->C = auglag->Ce;
230     }
231     if (!auglag->Y) {
232       auglag->Y = auglag->Ye;
233     }
234   }
235   /* initialize parameters */
236   if (auglag->type == TAO_ALMM_PHR) {
237     auglag->mu_fac = 10.0;
238     auglag->yi_min = 0.0;
239     auglag->ytol0 = 0.5;
240     auglag->gtol0 = tao->gatol;
241     if (tao->gatol_changed && tao->catol_changed) {
242       PetscCall(PetscInfo(tao,"TAOALMM with PHR: different gradient and constraint tolerances are not supported, setting catol = gatol\n"));
243       tao->catol = tao->gatol;
244     }
245   }
246   /* set the Lagrangian formulation type for the subsolver */
247   switch (auglag->type) {
248     case TAO_ALMM_CLASSIC:
249       auglag->sub_obj = TaoALMMComputeAugLagAndGradient_Private;
250       break;
251     case TAO_ALMM_PHR:
252       auglag->sub_obj = TaoALMMComputePHRLagAndGradient_Private;
253       break;
254     default:
255       break;
256   }
257   /* set up the subsolver */
258   PetscCall(TaoSetSolution(auglag->subsolver, auglag->P));
259   PetscCall(TaoSetObjective(auglag->subsolver, TaoALMMSubsolverObjective_Private, (void*)auglag));
260   PetscCall(TaoSetObjectiveAndGradient(auglag->subsolver, NULL, TaoALMMSubsolverObjectiveAndGradient_Private, (void*)auglag));
261   if (tao->bounded) {
262     /* make sure that the subsolver is a bound-constrained method */
263     PetscCall(PetscObjectTypeCompare((PetscObject)auglag->subsolver, TAOCG, &is_cg));
264     PetscCall(PetscObjectTypeCompare((PetscObject)auglag->subsolver, TAOLMVM, &is_lmvm));
265     if (is_cg) {
266       PetscCall(TaoSetType(auglag->subsolver, TAOBNCG));
267       PetscCall(PetscInfo(tao,"TAOCG detected for bound-constrained problem, switching to TAOBNCG instead."));
268     }
269     if (is_lmvm) {
270       PetscCall(TaoSetType(auglag->subsolver, TAOBQNLS));
271       PetscCall(PetscInfo(tao,"TAOLMVM detected for bound-constrained problem, switching to TAOBQNLS instead."));
272     }
273     /* create lower and upper bound clone vectors for subsolver */
274     if (!auglag->PL) {
275       PetscCall(VecDuplicate(auglag->P, &auglag->PL));
276     }
277     if (!auglag->PU) {
278       PetscCall(VecDuplicate(auglag->P, &auglag->PU));
279     }
280     if (tao->ineq_constrained) {
281       /* create lower and upper bounds for slack, set lower to 0 */
282       PetscCall(VecDuplicate(auglag->Ci, &SL));
283       PetscCall(VecSet(SL, 0.0));
284       PetscCall(VecDuplicate(auglag->Ci, &SU));
285       PetscCall(VecSet(SU, PETSC_INFINITY));
286       /* combine opt var bounds with slack bounds */
287       PetscCall(TaoALMMCombinePrimal_Private(tao, tao->XL, SL, auglag->PL));
288       PetscCall(TaoALMMCombinePrimal_Private(tao, tao->XU, SU, auglag->PU));
289       /* destroy work vectors */
290       PetscCall(VecDestroy(&SL));
291       PetscCall(VecDestroy(&SU));
292     } else {
293       /* no inequality constraints, just copy bounds into the subsolver */
294       PetscCall(VecCopy(tao->XL, auglag->PL));
295       PetscCall(VecCopy(tao->XU, auglag->PU));
296     }
297     PetscCall(TaoSetVariableBounds(auglag->subsolver, auglag->PL, auglag->PU));
298   }
299   PetscCall(TaoSetUp(auglag->subsolver));
300   auglag->G = auglag->subsolver->gradient;
301 
302   PetscFunctionReturn(0);
303 }
304 
305 static PetscErrorCode TaoDestroy_ALMM(Tao tao)
306 {
307   TAO_ALMM       *auglag = (TAO_ALMM*)tao->data;
308 
309   PetscFunctionBegin;
310   PetscCall(TaoDestroy(&auglag->subsolver));
311   if (tao->setupcalled) {
312     PetscCall(VecDestroy(&auglag->Xwork));              /* opt work */
313     if (tao->eq_constrained) {
314       PetscCall(VecDestroy(&auglag->Ye));               /* equality multipliers */
315       PetscCall(VecDestroy(&auglag->Cework));           /* equality work vector */
316     }
317     if (tao->ineq_constrained) {
318       PetscCall(VecDestroy(&auglag->Ps));               /* slack vars */
319       auglag->Parr[0] = NULL;                                     /* clear pointer to tao->solution, will be destroyed by TaoDestroy() shell */
320       PetscCall(PetscFree(auglag->Parr));               /* array of primal vectors */
321       PetscCall(VecDestroy(&auglag->LgradS));           /* slack grad */
322       PetscCall(VecDestroy(&auglag->Cizero));           /* zero vector for pointwise max */
323       PetscCall(VecDestroy(&auglag->Yi));               /* inequality multipliers */
324       PetscCall(VecDestroy(&auglag->Ciwork));           /* inequality work vector */
325       PetscCall(VecDestroy(&auglag->P));                /* combo primal */
326       PetscCall(ISDestroy(&auglag->Pis[0]));            /* index set for X inside P */
327       PetscCall(ISDestroy(&auglag->Pis[1]));            /* index set for S inside P */
328       PetscCall(PetscFree(auglag->Pis));                /* array of P index sets */
329       PetscCall(VecScatterDestroy(&auglag->Pscatter[0]));
330       PetscCall(VecScatterDestroy(&auglag->Pscatter[1]));
331       PetscCall(PetscFree(auglag->Pscatter));
332       if (tao->eq_constrained) {
333         PetscCall(VecDestroy(&auglag->Y));              /* combo multipliers */
334         PetscCall(PetscFree(auglag->Yarr));             /* array of dual vectors */
335         PetscCall(VecDestroy(&auglag->C));              /* combo constraints */
336         PetscCall(ISDestroy(&auglag->Yis[0]));          /* index set for Ye inside Y */
337         PetscCall(ISDestroy(&auglag->Yis[1]));          /* index set for Yi inside Y */
338         PetscCall(PetscFree(auglag->Yis));
339         PetscCall(VecScatterDestroy(&auglag->Yscatter[0]));
340         PetscCall(VecScatterDestroy(&auglag->Yscatter[1]));
341         PetscCall(PetscFree(auglag->Yscatter));
342       }
343     }
344     if (tao->bounded) {
345       PetscCall(VecDestroy(&auglag->PL));                /* lower bounds for subsolver */
346       PetscCall(VecDestroy(&auglag->PU));                /* upper bounds for subsolver */
347     }
348   }
349   PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoALMMGetType_C",NULL));
350   PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoALMMSetType_C",NULL));
351   PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoALMMGetSubsolver_C",NULL));
352   PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoALMMSetSubsolver_C",NULL));
353   PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoALMMGetMultipliers_C",NULL));
354   PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoALMMSetMultipliers_C",NULL));
355   PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoALMMGetPrimalIS_C",NULL));
356   PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoALMMGetDualIS_C",NULL));
357   PetscCall(PetscFree(tao->data));
358   PetscFunctionReturn(0);
359 }
360 
361 static PetscErrorCode TaoSetFromOptions_ALMM(PetscOptionItems *PetscOptionsObject,Tao tao)
362 {
363   TAO_ALMM       *auglag = (TAO_ALMM*)tao->data;
364   PetscInt       i;
365 
366   PetscFunctionBegin;
367   PetscOptionsHeadBegin(PetscOptionsObject,"Augmented Lagrangian multiplier method solves problems with general constraints by converting them into a sequence of unconstrained problems.");
368   PetscCall(PetscOptionsReal("-tao_almm_mu_init","initial penalty parameter","",auglag->mu0,&auglag->mu0,NULL));
369   PetscCall(PetscOptionsReal("-tao_almm_mu_factor","increase factor for the penalty parameter","",auglag->mu_fac,&auglag->mu_fac,NULL));
370   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));
371   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));
372   PetscCall(PetscOptionsReal("-tao_almm_mu_max","maximum safeguard for penalty parameter updates","",auglag->mu_max,&auglag->mu_max,NULL));
373   PetscCall(PetscOptionsReal("-tao_almm_ye_min","minimum safeguard for equality multiplier updates","",auglag->ye_min,&auglag->ye_min,NULL));
374   PetscCall(PetscOptionsReal("-tao_almm_ye_max","maximum safeguard for equality multipliers updates","",auglag->ye_max,&auglag->ye_max,NULL));
375   PetscCall(PetscOptionsReal("-tao_almm_yi_min","minimum safeguard for inequality multipliers updates","",auglag->yi_min,&auglag->yi_min,NULL));
376   PetscCall(PetscOptionsReal("-tao_almm_yi_max","maximum safeguard for inequality multipliers updates","",auglag->yi_max,&auglag->yi_max,NULL));
377   PetscCall(PetscOptionsEnum("-tao_almm_type","augmented Lagrangian formulation type for the subproblem","TaoALMMType",TaoALMMTypes,(PetscEnum)auglag->type,(PetscEnum*)&auglag->type,NULL));
378   PetscOptionsHeadEnd();
379   PetscCall(TaoSetOptionsPrefix(auglag->subsolver,((PetscObject)tao)->prefix));
380   PetscCall(TaoAppendOptionsPrefix(auglag->subsolver,"tao_almm_subsolver_"));
381   PetscCall(TaoSetFromOptions(auglag->subsolver));
382   for (i=0; i<tao->numbermonitors; i++) {
383     PetscCall(PetscObjectReference((PetscObject)tao->monitorcontext[i]));
384     PetscCall(TaoSetMonitor(auglag->subsolver, tao->monitor[i], tao->monitorcontext[i], tao->monitordestroy[i]));
385     if (tao->monitor[i] == TaoMonitorDefault || tao->monitor[i] == TaoDefaultCMonitor || tao->monitor[i] == TaoDefaultGMonitor || tao->monitor[i] == TaoDefaultSMonitor) {
386       auglag->info = PETSC_TRUE;
387     }
388   }
389   PetscFunctionReturn(0);
390 }
391 
392 /* -------------------------------------------------------- */
393 
394 /*MC
395   TaoALMM - Augmented Lagrangian multiplier method for solving nonlinear optimization problems with general constraints.
396 
397   Options Database Keys:
398 + -tao_almm_mu_init <real>       - initial penalty parameter (default: 10.)
399 . -tao_almm_mu_factor <real>     - increase factor for the penalty parameter (default: 100.)
400 . -tao_almm_mu_max <real>        - maximum safeguard for penalty parameter updates (default: 1.e20)
401 . -tao_almm_mu_power_good <real> - exponential for penalty parameter when multiplier update is accepted (default: 0.9)
402 . -tao_almm_mu_power_bad <real>  - exponential for penalty parameter when multiplier update is rejected (default: 0.1)
403 . -tao_almm_ye_min <real>        - minimum safeguard for equality multiplier updates (default: -1.e20)
404 . -tao_almm_ye_max <real>        - maximum safeguard for equality multiplier updates (default: 1.e20)
405 . -tao_almm_yi_min <real>        - minimum safeguard for inequality multiplier updates (default: -1.e20)
406 . -tao_almm_yi_max <real>        - maximum safeguard for inequality multiplier updates (default: 1.e20)
407 - -tao_almm_type <classic,phr>   - change formulation of the augmented Lagrangian merit function for the subproblem (default: classic)
408 
409   Level: beginner
410 
411   Notes:
412   This method converts a constrained problem into a sequence of unconstrained problems via the augmented
413   Lagrangian merit function. Bound constraints are pushed down to the subproblem without any modifications.
414 
415   Two formulations are offered for the subproblem: canonical Hestenes-Powell augmented Lagrangian with slack
416   variables for inequality constraints, and a slack-less Powell-Hestenes-Rockafellar (PHR) formulation utilizing a
417   pointwise max() penalty on inequality constraints. The canonical augmented Lagrangian formulation typically
418   converges faster for most problems. However, PHR may be desirable for problems featuring a large number
419   of inequality constraints because it avoids inflating the size of the subproblem with slack variables.
420 
421   The subproblem is solved using a nested first-order TAO solver. The user can retrieve a pointer to
422   the subsolver via TaoALMMGetSubsolver() or pass command line arguments to it using the
423   "-tao_almm_subsolver_" prefix. Currently, TaoALMM does not support second-order methods for the
424   subproblem. It is also highly recommended that the subsolver chosen by the user utilize a trust-region
425   strategy for globalization (default: TAOBQNKTR) especially if the outer problem features bound constraints.
426 
427 .vb
428   while unconverged
429     solve argmin_x L(x) s.t. l <= x <= u
430     if ||c|| <= y_tol
431       if ||c|| <= c_tol && ||Lgrad|| <= g_tol:
432         problem converged, return solution
433       else
434         constraints sufficiently improved
435         update multipliers and tighten tolerances
436       endif
437     else
438       constraints did not improve
439       update penalty and loosen tolerances
440     endif
441   endwhile
442 .ve
443 
444 .seealso: `TaoALMMGetType()`, `TaoALMMSetType()`, `TaoALMMSetSubsolver()`, `TaoALMMGetSubsolver()`,
445           `TaoALMMGetMultipliers()`, `TaoALMMSetMultipliers()`, `TaoALMMGetPrimalIS()`, `TaoALMMGetDualIS()`
446 M*/
447 PETSC_EXTERN PetscErrorCode TaoCreate_ALMM(Tao tao)
448 {
449   TAO_ALMM       *auglag;
450 
451   PetscFunctionBegin;
452   PetscCall(PetscNewLog(tao, &auglag));
453 
454   tao->ops->destroy        = TaoDestroy_ALMM;
455   tao->ops->setup          = TaoSetUp_ALMM;
456   tao->ops->setfromoptions = TaoSetFromOptions_ALMM;
457   tao->ops->view           = TaoView_ALMM;
458   tao->ops->solve          = TaoSolve_ALMM;
459 
460   tao->gatol = 1.e-5;
461   tao->grtol = 0.0;
462   tao->gttol = 0.0;
463   tao->catol = 1.e-5;
464   tao->crtol = 0.0;
465 
466   tao->data           = (void*)auglag;
467   auglag->parent      = tao;
468   auglag->mu0         = 10.0;
469   auglag->mu          = auglag->mu0;
470   auglag->mu_fac      = 10.0;
471   auglag->mu_max      = PETSC_INFINITY;
472   auglag->mu_pow_good = 0.9;
473   auglag->mu_pow_bad  = 0.1;
474   auglag->ye_min      = PETSC_NINFINITY;
475   auglag->ye_max      = PETSC_INFINITY;
476   auglag->yi_min      = PETSC_NINFINITY;
477   auglag->yi_max      = PETSC_INFINITY;
478   auglag->ytol0       = 0.1/PetscPowReal(auglag->mu0, auglag->mu_pow_bad);
479   auglag->ytol        = auglag->ytol0;
480   auglag->gtol0       = 1.0/auglag->mu0;
481   auglag->gtol        = auglag->gtol0;
482 
483   auglag->sub_obj     = TaoALMMComputeAugLagAndGradient_Private;
484   auglag->type        = TAO_ALMM_CLASSIC;
485   auglag->info        = PETSC_FALSE;
486 
487   PetscCall(TaoCreate(PetscObjectComm((PetscObject)tao),&auglag->subsolver));
488   PetscCall(TaoSetType(auglag->subsolver, TAOBQNKTR));
489   PetscCall(TaoSetTolerances(auglag->subsolver, auglag->gtol, 0.0, 0.0));
490   PetscCall(TaoSetMaximumIterations(auglag->subsolver, 1000));
491   PetscCall(TaoSetMaximumFunctionEvaluations(auglag->subsolver, 10000));
492   PetscCall(TaoSetFunctionLowerBound(auglag->subsolver, PETSC_NINFINITY));
493   PetscCall(PetscObjectIncrementTabLevel((PetscObject)auglag->subsolver,(PetscObject)tao,1));
494 
495   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetType_C", TaoALMMGetType_Private));
496   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMSetType_C", TaoALMMSetType_Private));
497   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetSubsolver_C", TaoALMMGetSubsolver_Private));
498   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMSetSubsolver_C", TaoALMMSetSubsolver_Private));
499   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetMultipliers_C", TaoALMMGetMultipliers_Private));
500   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMSetMultipliers_C", TaoALMMSetMultipliers_Private));
501   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetPrimalIS_C", TaoALMMGetPrimalIS_Private));
502   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetDualIS_C", TaoALMMGetDualIS_Private));
503   PetscFunctionReturn(0);
504 }
505 
506 static PetscErrorCode TaoALMMCombinePrimal_Private(Tao tao, Vec X, Vec S, Vec P)
507 {
508   TAO_ALMM     *auglag = (TAO_ALMM*)tao->data;
509 
510   PetscFunctionBegin;
511   if (tao->ineq_constrained) {
512     PetscCall(VecScatterBegin(auglag->Pscatter[0], X, P, INSERT_VALUES, SCATTER_REVERSE));
513     PetscCall(VecScatterEnd(auglag->Pscatter[0], X, P, INSERT_VALUES, SCATTER_REVERSE));
514     PetscCall(VecScatterBegin(auglag->Pscatter[1], S, P, INSERT_VALUES, SCATTER_REVERSE));
515     PetscCall(VecScatterEnd(auglag->Pscatter[1], S, P, INSERT_VALUES, SCATTER_REVERSE));
516   } else {
517     PetscCall(VecCopy(X, P));
518   }
519   PetscFunctionReturn(0);
520 }
521 
522 static PetscErrorCode TaoALMMCombineDual_Private(Tao tao, Vec EQ, Vec IN, Vec Y)
523 {
524   TAO_ALMM     *auglag = (TAO_ALMM*)tao->data;
525 
526   PetscFunctionBegin;
527   if (tao->eq_constrained) {
528     if (tao->ineq_constrained) {
529       PetscCall(VecScatterBegin(auglag->Yscatter[0], EQ, Y, INSERT_VALUES, SCATTER_REVERSE));
530       PetscCall(VecScatterEnd(auglag->Yscatter[0], EQ, Y, INSERT_VALUES, SCATTER_REVERSE));
531       PetscCall(VecScatterBegin(auglag->Yscatter[1], IN, Y, INSERT_VALUES, SCATTER_REVERSE));
532       PetscCall(VecScatterEnd(auglag->Yscatter[1], IN, Y, INSERT_VALUES, SCATTER_REVERSE));
533     } else {
534       PetscCall(VecCopy(EQ, Y));
535     }
536   } else {
537     PetscCall(VecCopy(IN, Y));
538   }
539   PetscFunctionReturn(0);
540 }
541 
542 static PetscErrorCode TaoALMMSplitPrimal_Private(Tao tao, Vec P, Vec X, Vec S)
543 {
544   TAO_ALMM     *auglag = (TAO_ALMM*)tao->data;
545 
546   PetscFunctionBegin;
547   if (tao->ineq_constrained) {
548     PetscCall(VecScatterBegin(auglag->Pscatter[0], P, X, INSERT_VALUES, SCATTER_FORWARD));
549     PetscCall(VecScatterEnd(auglag->Pscatter[0], P, X, INSERT_VALUES, SCATTER_FORWARD));
550     PetscCall(VecScatterBegin(auglag->Pscatter[1], P, S, INSERT_VALUES, SCATTER_FORWARD));
551     PetscCall(VecScatterEnd(auglag->Pscatter[1], P, S, INSERT_VALUES, SCATTER_FORWARD));
552   } else {
553     PetscCall(VecCopy(P, X));
554   }
555   PetscFunctionReturn(0);
556 }
557 
558 /* this assumes that the latest constraints are stored in Ce and Ci, and also combined in C */
559 static PetscErrorCode TaoALMMComputeOptimalityNorms_Private(Tao tao)
560 {
561   TAO_ALMM     *auglag = (TAO_ALMM*)tao->data;
562 
563   PetscFunctionBegin;
564   /* if bounded, project the gradient */
565   if (tao->bounded) PetscCall(VecBoundGradientProjection(auglag->LgradX, auglag->Px, tao->XL, tao->XU, auglag->LgradX));
566   if (auglag->type == TAO_ALMM_PHR) {
567     PetscCall(VecNorm(auglag->LgradX, NORM_INFINITY, &auglag->gnorm));
568     auglag->cenorm = 0.0;
569     if (tao->eq_constrained) {
570       PetscCall(VecNorm(auglag->Ce, NORM_INFINITY, &auglag->cenorm));
571     }
572     auglag->cinorm = 0.0;
573     if (tao->ineq_constrained) {
574       PetscCall(VecCopy(auglag->Yi, auglag->Ciwork));
575       PetscCall(VecScale(auglag->Ciwork, -1.0/auglag->mu));
576       PetscCall(VecPointwiseMax(auglag->Ciwork, auglag->Ci, auglag->Ciwork));
577       PetscCall(VecNorm(auglag->Ciwork, NORM_INFINITY, &auglag->cinorm));
578     }
579     auglag->cnorm_old = auglag->cnorm;
580     auglag->cnorm = PetscMax(auglag->cenorm, auglag->cinorm);
581     auglag->ytol = auglag->ytol0 * auglag->cnorm_old;
582   } else {
583     PetscCall(VecNorm(auglag->LgradX, NORM_2, &auglag->gnorm));
584     PetscCall(VecNorm(auglag->C, NORM_2, &auglag->cnorm));
585   }
586   PetscFunctionReturn(0);
587 }
588 
589 static PetscErrorCode TaoALMMEvaluateIterate_Private(Tao tao, Vec P)
590 {
591   TAO_ALMM       *auglag = (TAO_ALMM*)tao->data;
592 
593   PetscFunctionBegin;
594   /* split solution into primal and slack components */
595   PetscCall(TaoALMMSplitPrimal_Private(tao, auglag->P, auglag->Px, auglag->Ps));
596 
597   /* compute f, df/dx and the constraints */
598   PetscCall(TaoComputeObjectiveAndGradient(tao, auglag->Px, &auglag->fval, auglag->LgradX));
599   if (tao->eq_constrained) {
600     PetscCall(TaoComputeEqualityConstraints(tao, auglag->Px, auglag->Ce));
601     PetscCall(TaoComputeJacobianEquality(tao, auglag->Px, auglag->Ae, auglag->Ae));
602   }
603   if (tao->ineq_constrained) {
604     PetscCall(TaoComputeInequalityConstraints(tao, auglag->Px, auglag->Ci));
605     PetscCall(TaoComputeJacobianInequality(tao, auglag->Px, auglag->Ai, auglag->Ai));
606     switch (auglag->type) {
607       case TAO_ALMM_CLASSIC:
608         /* classic formulation converts inequality to equality constraints via slack variables */
609         PetscCall(VecAXPY(auglag->Ci, -1.0, auglag->Ps));
610         break;
611       case TAO_ALMM_PHR:
612         /* PHR is based on Ci <= 0 while TAO defines Ci >= 0 so we hit it with a negative sign */
613         PetscCall(VecScale(auglag->Ci, -1.0));
614         PetscCall(MatScale(auglag->Ai, -1.0));
615         break;
616       default:
617         break;
618     }
619   }
620   /* combine constraints into one vector */
621   PetscCall(TaoALMMCombineDual_Private(tao, auglag->Ce, auglag->Ci, auglag->C));
622   PetscFunctionReturn(0);
623 }
624 
625 /*
626 Lphr = f + 0.5*mu*[ (Ce + Ye/mu)^T (Ce + Ye/mu) + pmin(0, Ci + Yi/mu)^T pmin(0, Ci + Yi/mu)]
627 
628 dLphr/dX = dF/dX + mu*[ (Ce + Ye/mu)^T Ae + pmin(0, Ci + Yi/mu)^T Ai]
629 
630 dLphr/dS = 0
631 */
632 static PetscErrorCode TaoALMMComputePHRLagAndGradient_Private(Tao tao)
633 {
634   TAO_ALMM       *auglag = (TAO_ALMM*)tao->data;
635   PetscReal      eq_norm=0.0, ineq_norm=0.0;
636 
637   PetscFunctionBegin;
638   PetscCall(TaoALMMEvaluateIterate_Private(tao, auglag->P));
639   if (tao->eq_constrained) {
640     /* Ce_work = mu*(Ce + Ye/mu) */
641     PetscCall(VecWAXPY(auglag->Cework, 1.0/auglag->mu, auglag->Ye, auglag->Ce));
642     PetscCall(VecDot(auglag->Cework, auglag->Cework, &eq_norm)); /* contribution to scalar Lagrangian */
643     PetscCall(VecScale(auglag->Cework, auglag->mu));
644     /* dL/dX += mu*(Ce + Ye/mu)^T Ae */
645     PetscCall(MatMultTransposeAdd(auglag->Ae, auglag->Cework, auglag->LgradX, auglag->LgradX));
646   }
647   if (tao->ineq_constrained) {
648     /* Ci_work = mu * pmax(0, Ci + Yi/mu) where pmax() is pointwise max() */
649     PetscCall(VecWAXPY(auglag->Ciwork, 1.0/auglag->mu, auglag->Yi, auglag->Ci));
650     PetscCall(VecPointwiseMax(auglag->Ciwork, auglag->Cizero, auglag->Ciwork));
651     PetscCall(VecDot(auglag->Ciwork, auglag->Ciwork, &ineq_norm)); /* contribution to scalar Lagrangian */
652     /* dL/dX += mu * pmax(0, Ci + Yi/mu)^T Ai */
653     PetscCall(VecScale(auglag->Ciwork, auglag->mu));
654     PetscCall(MatMultTransposeAdd(auglag->Ai, auglag->Ciwork, auglag->LgradX, auglag->LgradX));
655     /* dL/dS = 0 because there are no slacks in PHR */
656     PetscCall(VecZeroEntries(auglag->LgradS));
657   }
658   /* combine gradient together */
659   PetscCall(TaoALMMCombinePrimal_Private(tao, auglag->LgradX, auglag->LgradS, auglag->G));
660   /* compute L = f + 0.5 * mu * [(Ce + Ye/mu)^T (Ce + Ye/mu) + pmax(0, Ci + Yi/mu)^T pmax(0, Ci + Yi/mu)] */
661   auglag->Lval = auglag->fval + 0.5*auglag->mu*(eq_norm + ineq_norm);
662   PetscFunctionReturn(0);
663 }
664 
665 /*
666 Lc = F + Ye^TCe + Yi^T(Ci - S) + 0.5*mu*[Ce^TCe + (Ci - S)^T(Ci - S)]
667 
668 dLc/dX = dF/dX + Ye^TAe + Yi^TAi + 0.5*mu*[Ce^TAe + (Ci - S)^TAi]
669 
670 dLc/dS = -[Yi + mu*(Ci - S)]
671 */
672 static PetscErrorCode TaoALMMComputeAugLagAndGradient_Private(Tao tao)
673 {
674   TAO_ALMM       *auglag = (TAO_ALMM*)tao->data;
675   PetscReal      yeTce=0.0, yiTcims=0.0, ceTce=0.0, cimsTcims=0.0;
676 
677   PetscFunctionBegin;
678   PetscCall(TaoALMMEvaluateIterate_Private(tao, auglag->P));
679   if (tao->eq_constrained) {
680     /* compute scalar contributions */
681     PetscCall(VecDot(auglag->Ye, auglag->Ce, &yeTce));
682     PetscCall(VecDot(auglag->Ce, auglag->Ce, &ceTce));
683     /* dL/dX += ye^T Ae */
684     PetscCall(MatMultTransposeAdd(auglag->Ae, auglag->Ye, auglag->LgradX, auglag->LgradX));
685     /* dL/dX += mu * ce^T Ae */
686     PetscCall(MatMultTranspose(auglag->Ae, auglag->Ce, auglag->Xwork));
687     PetscCall(VecAXPY(auglag->LgradX, auglag->mu, auglag->Xwork));
688   }
689   if (tao->ineq_constrained) {
690     /* compute scalar contributions */
691     PetscCall(VecDot(auglag->Yi, auglag->Ci, &yiTcims));
692     PetscCall(VecDot(auglag->Ci, auglag->Ci, &cimsTcims));
693     /* dL/dX += yi^T Ai */
694     PetscCall(MatMultTransposeAdd(auglag->Ai, auglag->Yi, auglag->LgradX, auglag->LgradX));
695     /* dL/dX += mu * (ci - s)^T Ai */
696     PetscCall(MatMultTranspose(auglag->Ai, auglag->Ci, auglag->Xwork));
697     PetscCall(VecAXPY(auglag->LgradX, auglag->mu, auglag->Xwork));
698     /* dL/dS = -[yi + mu*(ci - s)] */
699     PetscCall(VecWAXPY(auglag->LgradS, auglag->mu, auglag->Ci, auglag->Yi));
700     PetscCall(VecScale(auglag->LgradS, -1.0));
701   }
702   /* combine gradient together */
703   PetscCall(TaoALMMCombinePrimal_Private(tao, auglag->LgradX, auglag->LgradS, auglag->G));
704   /* compute L = f + ye^T ce + yi^T (ci - s) + 0.5*mu*||ce||^2 + 0.5*mu*||ci - s||^2 */
705   auglag->Lval = auglag->fval + yeTce + yiTcims + 0.5*auglag->mu*(ceTce + cimsTcims);
706   PetscFunctionReturn(0);
707 }
708 
709 PetscErrorCode TaoALMMSubsolverObjective_Private(Tao tao, Vec P, PetscReal *Lval, void *ctx)
710 {
711   TAO_ALMM       *auglag = (TAO_ALMM*)ctx;
712 
713   PetscFunctionBegin;
714   PetscCall(VecCopy(P, auglag->P));
715   PetscCall((*auglag->sub_obj)(auglag->parent));
716   *Lval = auglag->Lval;
717   PetscFunctionReturn(0);
718 }
719 
720 PetscErrorCode TaoALMMSubsolverObjectiveAndGradient_Private(Tao tao, Vec P, PetscReal *Lval, Vec G, void *ctx)
721 {
722   TAO_ALMM       *auglag = (TAO_ALMM*)ctx;
723 
724   PetscFunctionBegin;
725   PetscCall(VecCopy(P, auglag->P));
726   PetscCall((*auglag->sub_obj)(auglag->parent));
727   PetscCall(VecCopy(auglag->G, G));
728   *Lval = auglag->Lval;
729   PetscFunctionReturn(0);
730 }
731