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