xref: /petsc/src/tao/constrained/impls/almm/almm.c (revision 4dfa11a44d5adf2389f1d3acbc8f3c1116dc6c3a)
1661095bbSAlp Dener #include <../src/tao/constrained/impls/almm/almm.h> /*I "petsctao.h" I*/ /*I "petscvec.h" I*/
2661095bbSAlp Dener #include <petsctao.h>
3661095bbSAlp Dener #include <petsc/private/petscimpl.h>
4661095bbSAlp Dener #include <petsc/private/vecimpl.h>
5661095bbSAlp Dener 
6661095bbSAlp Dener static PetscErrorCode TaoALMMCombinePrimal_Private(Tao, Vec, Vec, Vec);
7661095bbSAlp Dener static PetscErrorCode TaoALMMCombineDual_Private(Tao, Vec, Vec, Vec);
8661095bbSAlp Dener static PetscErrorCode TaoALMMSplitPrimal_Private(Tao, Vec, Vec, Vec);
9661095bbSAlp Dener static PetscErrorCode TaoALMMComputeOptimalityNorms_Private(Tao);
10661095bbSAlp Dener static PetscErrorCode TaoALMMComputeAugLagAndGradient_Private(Tao);
11661095bbSAlp Dener static PetscErrorCode TaoALMMComputePHRLagAndGradient_Private(Tao);
12661095bbSAlp Dener 
139371c9d4SSatish Balay static PetscErrorCode TaoSolve_ALMM(Tao tao) {
14661095bbSAlp Dener   TAO_ALMM          *auglag = (TAO_ALMM *)tao->data;
15661095bbSAlp Dener   TaoConvergedReason reason;
16661095bbSAlp Dener   PetscReal          updated;
17661095bbSAlp Dener 
18661095bbSAlp Dener   PetscFunctionBegin;
19661095bbSAlp Dener   /* reset initial multiplier/slack guess */
20661095bbSAlp Dener   if (!tao->recycle) {
21661095bbSAlp Dener     if (tao->ineq_constrained) {
229566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(auglag->Ps));
239566063dSJacob Faibussowitsch       PetscCall(TaoALMMCombinePrimal_Private(tao, auglag->Px, auglag->Ps, auglag->P));
249566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(auglag->Yi));
25661095bbSAlp Dener     }
261baa6e33SBarry Smith     if (tao->eq_constrained) PetscCall(VecZeroEntries(auglag->Ye));
27661095bbSAlp Dener   }
28661095bbSAlp Dener 
29661095bbSAlp Dener   /* compute initial nonlinear Lagrangian and its derivatives */
309566063dSJacob Faibussowitsch   PetscCall((*auglag->sub_obj)(tao));
319566063dSJacob Faibussowitsch   PetscCall(TaoALMMComputeOptimalityNorms_Private(tao));
32661095bbSAlp Dener   /* print initial step and check convergence */
339566063dSJacob Faibussowitsch   PetscCall(PetscInfo(tao, "Solving with %s formulation\n", TaoALMMTypes[auglag->type]));
349566063dSJacob Faibussowitsch   PetscCall(TaoLogConvergenceHistory(tao, auglag->Lval, auglag->gnorm, auglag->cnorm, tao->ksp_its));
359566063dSJacob Faibussowitsch   PetscCall(TaoMonitor(tao, tao->niter, auglag->fval, auglag->gnorm, auglag->cnorm, 0.0));
36dbbe0bcdSBarry Smith   PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
37661095bbSAlp Dener   /* set initial penalty factor and inner solver tolerance */
38661095bbSAlp Dener   switch (auglag->type) {
399371c9d4SSatish Balay   case TAO_ALMM_CLASSIC: auglag->mu = auglag->mu0; break;
407721a36fSStefano Zampini   case TAO_ALMM_PHR:
41661095bbSAlp Dener     auglag->cenorm = 0.0;
4248a46eb9SPierre Jolivet     if (tao->eq_constrained) PetscCall(VecDot(auglag->Ce, auglag->Ce, &auglag->cenorm));
43661095bbSAlp Dener     auglag->cinorm = 0.0;
44661095bbSAlp Dener     if (tao->ineq_constrained) {
459566063dSJacob Faibussowitsch       PetscCall(VecCopy(auglag->Ci, auglag->Ciwork));
469566063dSJacob Faibussowitsch       PetscCall(VecScale(auglag->Ciwork, -1.0));
479566063dSJacob Faibussowitsch       PetscCall(VecPointwiseMax(auglag->Ciwork, auglag->Cizero, auglag->Ciwork));
489566063dSJacob Faibussowitsch       PetscCall(VecDot(auglag->Ciwork, auglag->Ciwork, &auglag->cinorm));
49661095bbSAlp Dener     }
50661095bbSAlp Dener     /* determine initial penalty factor based on the balance of constraint violation and objective function value */
51661095bbSAlp Dener     auglag->mu = PetscMax(1.e-6, PetscMin(10.0, 2.0 * PetscAbsReal(auglag->fval) / (auglag->cenorm + auglag->cinorm)));
52661095bbSAlp Dener     break;
539371c9d4SSatish Balay   default: break;
54661095bbSAlp Dener   }
55661095bbSAlp Dener   auglag->gtol = auglag->gtol0;
5663a3b9bcSJacob Faibussowitsch   PetscCall(PetscInfo(tao, "Initial penalty: %.2g\n", (double)auglag->mu));
57661095bbSAlp Dener 
58661095bbSAlp Dener   /* start aug-lag outer loop */
59661095bbSAlp Dener   while (tao->reason == TAO_CONTINUE_ITERATING) {
60661095bbSAlp Dener     ++tao->niter;
61661095bbSAlp Dener     /* update subsolver tolerance */
6263a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(tao, "Subsolver tolerance: ||G|| <= %e\n", (double)auglag->gtol));
639566063dSJacob Faibussowitsch     PetscCall(TaoSetTolerances(auglag->subsolver, auglag->gtol, 0.0, 0.0));
64661095bbSAlp Dener     /* solve the bound-constrained or unconstrained subproblem */
659566063dSJacob Faibussowitsch     PetscCall(TaoSolve(auglag->subsolver));
669566063dSJacob Faibussowitsch     PetscCall(TaoGetConvergedReason(auglag->subsolver, &reason));
67661095bbSAlp Dener     tao->ksp_its += auglag->subsolver->ksp_its;
6848a46eb9SPierre Jolivet     if (reason != TAO_CONVERGED_GATOL) PetscCall(PetscInfo(tao, "Subsolver failed to converge, reason: %s\n", TaoConvergedReasons[reason]));
69661095bbSAlp Dener     /* evaluate solution and test convergence */
709566063dSJacob Faibussowitsch     PetscCall((*auglag->sub_obj)(tao));
719566063dSJacob Faibussowitsch     PetscCall(TaoALMMComputeOptimalityNorms_Private(tao));
72661095bbSAlp Dener     /* decide whether to update multipliers or not */
73661095bbSAlp Dener     updated = 0.0;
74661095bbSAlp Dener     if (auglag->cnorm <= auglag->ytol) {
7563a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(tao, "Multipliers updated: ||C|| <= %e\n", (double)auglag->ytol));
76661095bbSAlp Dener       /* constraints are good, update multipliers and convergence tolerances */
77661095bbSAlp Dener       if (tao->eq_constrained) {
789566063dSJacob Faibussowitsch         PetscCall(VecAXPY(auglag->Ye, auglag->mu, auglag->Ce));
799566063dSJacob Faibussowitsch         PetscCall(VecSet(auglag->Cework, auglag->ye_max));
809566063dSJacob Faibussowitsch         PetscCall(VecPointwiseMin(auglag->Ye, auglag->Cework, auglag->Ye));
819566063dSJacob Faibussowitsch         PetscCall(VecSet(auglag->Cework, auglag->ye_min));
829566063dSJacob Faibussowitsch         PetscCall(VecPointwiseMax(auglag->Ye, auglag->Cework, auglag->Ye));
83661095bbSAlp Dener       }
84661095bbSAlp Dener       if (tao->ineq_constrained) {
859566063dSJacob Faibussowitsch         PetscCall(VecAXPY(auglag->Yi, auglag->mu, auglag->Ci));
869566063dSJacob Faibussowitsch         PetscCall(VecSet(auglag->Ciwork, auglag->yi_max));
879566063dSJacob Faibussowitsch         PetscCall(VecPointwiseMin(auglag->Yi, auglag->Ciwork, auglag->Yi));
889566063dSJacob Faibussowitsch         PetscCall(VecSet(auglag->Ciwork, auglag->yi_min));
899566063dSJacob Faibussowitsch         PetscCall(VecPointwiseMax(auglag->Yi, auglag->Ciwork, auglag->Yi));
90661095bbSAlp Dener       }
91661095bbSAlp Dener       /* tolerances are updated only for non-PHR methods */
92661095bbSAlp Dener       if (auglag->type != TAO_ALMM_PHR) {
93661095bbSAlp Dener         auglag->ytol = PetscMax(tao->catol, auglag->ytol / PetscPowReal(auglag->mu, auglag->mu_pow_good));
94661095bbSAlp Dener         auglag->gtol = PetscMax(tao->gatol, auglag->gtol / auglag->mu);
95661095bbSAlp Dener       }
96661095bbSAlp Dener       updated = 1.0;
97661095bbSAlp Dener     } else {
98661095bbSAlp Dener       /* constraints are bad, update penalty factor */
99661095bbSAlp Dener       auglag->mu = PetscMin(auglag->mu_max, auglag->mu_fac * auglag->mu);
100661095bbSAlp Dener       /* tolerances are reset only for non-PHR methods */
101661095bbSAlp Dener       if (auglag->type != TAO_ALMM_PHR) {
102661095bbSAlp Dener         auglag->ytol = PetscMax(tao->catol, 0.1 / PetscPowReal(auglag->mu, auglag->mu_pow_bad));
103661095bbSAlp Dener         auglag->gtol = PetscMax(tao->gatol, 1.0 / auglag->mu);
104661095bbSAlp Dener       }
10563a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(tao, "Penalty increased: mu = %.2g\n", (double)auglag->mu));
106661095bbSAlp Dener     }
1079566063dSJacob Faibussowitsch     PetscCall(TaoLogConvergenceHistory(tao, auglag->fval, auglag->gnorm, auglag->cnorm, tao->ksp_its));
1089566063dSJacob Faibussowitsch     PetscCall(TaoMonitor(tao, tao->niter, auglag->fval, auglag->gnorm, auglag->cnorm, updated));
109dbbe0bcdSBarry Smith     PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
110661095bbSAlp Dener   }
111661095bbSAlp Dener 
112661095bbSAlp Dener   PetscFunctionReturn(0);
113661095bbSAlp Dener }
114661095bbSAlp Dener 
1159371c9d4SSatish Balay static PetscErrorCode TaoView_ALMM(Tao tao, PetscViewer viewer) {
116661095bbSAlp Dener   TAO_ALMM *auglag = (TAO_ALMM *)tao->data;
117661095bbSAlp Dener   PetscBool isascii;
118661095bbSAlp Dener 
119661095bbSAlp Dener   PetscFunctionBegin;
1209566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPushTab(viewer));
1219566063dSJacob Faibussowitsch   PetscCall(TaoView(auglag->subsolver, viewer));
1229566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPopTab(viewer));
1239566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
124661095bbSAlp Dener   if (isascii) {
1259566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
1269566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "ALMM Formulation Type: %s\n", TaoALMMTypes[auglag->type]));
1279566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
128661095bbSAlp Dener   }
129661095bbSAlp Dener   PetscFunctionReturn(0);
130661095bbSAlp Dener }
131661095bbSAlp Dener 
1329371c9d4SSatish Balay static PetscErrorCode TaoSetUp_ALMM(Tao tao) {
133661095bbSAlp Dener   TAO_ALMM *auglag = (TAO_ALMM *)tao->data;
134661095bbSAlp Dener   VecType   vec_type;
135661095bbSAlp Dener   Vec       SL, SU;
136661095bbSAlp Dener   PetscBool is_cg = PETSC_FALSE, is_lmvm = PETSC_FALSE;
137661095bbSAlp Dener 
138661095bbSAlp Dener   PetscFunctionBegin;
1393c859ba3SBarry Smith   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.");
1403c859ba3SBarry Smith   PetscCheck(tao->eq_constrained || tao->ineq_constrained, PetscObjectComm((PetscObject)tao), PETSC_ERR_ORDER, "Equality and/or inequality constraints must be defined before solver setup.");
1419566063dSJacob Faibussowitsch   PetscCall(TaoComputeVariableBounds(tao));
142661095bbSAlp Dener   /* alias base vectors and create extras */
1439566063dSJacob Faibussowitsch   PetscCall(VecGetType(tao->solution, &vec_type));
144661095bbSAlp Dener   auglag->Px = tao->solution;
145661095bbSAlp Dener   if (!tao->gradient) { /* base gradient */
1469566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(tao->solution, &tao->gradient));
147661095bbSAlp Dener   }
148661095bbSAlp Dener   auglag->LgradX = tao->gradient;
149661095bbSAlp Dener   if (!auglag->Xwork) { /* opt var work vector */
1509566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(tao->solution, &auglag->Xwork));
151661095bbSAlp Dener   }
152661095bbSAlp Dener   if (tao->eq_constrained) {
153661095bbSAlp Dener     auglag->Ce = tao->constraints_equality;
154661095bbSAlp Dener     auglag->Ae = tao->jacobian_equality;
155661095bbSAlp Dener     if (!auglag->Ye) { /* equality multipliers */
1569566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(auglag->Ce, &auglag->Ye));
157661095bbSAlp Dener     }
15848a46eb9SPierre Jolivet     if (!auglag->Cework) PetscCall(VecDuplicate(auglag->Ce, &auglag->Cework));
159661095bbSAlp Dener   }
160661095bbSAlp Dener   if (tao->ineq_constrained) {
161661095bbSAlp Dener     auglag->Ci = tao->constraints_inequality;
162661095bbSAlp Dener     auglag->Ai = tao->jacobian_inequality;
163661095bbSAlp Dener     if (!auglag->Yi) { /* inequality multipliers */
1649566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(auglag->Ci, &auglag->Yi));
165661095bbSAlp Dener     }
16648a46eb9SPierre Jolivet     if (!auglag->Ciwork) PetscCall(VecDuplicate(auglag->Ci, &auglag->Ciwork));
167661095bbSAlp Dener     if (!auglag->Cizero) {
1689566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(auglag->Ci, &auglag->Cizero));
1699566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(auglag->Cizero));
170661095bbSAlp Dener     }
171661095bbSAlp Dener     if (!auglag->Ps) { /* slack vars */
1729566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(auglag->Ci, &auglag->Ps));
173661095bbSAlp Dener     }
174661095bbSAlp Dener     if (!auglag->LgradS) { /* slack component of Lagrangian gradient */
1759566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(auglag->Ci, &auglag->LgradS));
176661095bbSAlp Dener     }
177661095bbSAlp Dener     /* create vector for combined primal space and the associated communication objects */
178661095bbSAlp Dener     if (!auglag->P) {
1799566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(2, &auglag->Parr));
1809371c9d4SSatish Balay       auglag->Parr[0] = auglag->Px;
1819371c9d4SSatish Balay       auglag->Parr[1] = auglag->Ps;
1829566063dSJacob Faibussowitsch       PetscCall(VecConcatenate(2, auglag->Parr, &auglag->P, &auglag->Pis));
1839566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(2, &auglag->Pscatter));
1849566063dSJacob Faibussowitsch       PetscCall(VecScatterCreate(auglag->P, auglag->Pis[0], auglag->Px, NULL, &auglag->Pscatter[0]));
1859566063dSJacob Faibussowitsch       PetscCall(VecScatterCreate(auglag->P, auglag->Pis[1], auglag->Ps, NULL, &auglag->Pscatter[1]));
186661095bbSAlp Dener     }
187661095bbSAlp Dener     if (tao->eq_constrained) {
188661095bbSAlp Dener       /* create vector for combined dual space and the associated communication objects */
189661095bbSAlp Dener       if (!auglag->Y) {
1909566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(2, &auglag->Yarr));
1919371c9d4SSatish Balay         auglag->Yarr[0] = auglag->Ye;
1929371c9d4SSatish Balay         auglag->Yarr[1] = auglag->Yi;
1939566063dSJacob Faibussowitsch         PetscCall(VecConcatenate(2, auglag->Yarr, &auglag->Y, &auglag->Yis));
1949566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(2, &auglag->Yscatter));
1959566063dSJacob Faibussowitsch         PetscCall(VecScatterCreate(auglag->Y, auglag->Yis[0], auglag->Ye, NULL, &auglag->Yscatter[0]));
1969566063dSJacob Faibussowitsch         PetscCall(VecScatterCreate(auglag->Y, auglag->Yis[1], auglag->Yi, NULL, &auglag->Yscatter[1]));
197661095bbSAlp Dener       }
19848a46eb9SPierre Jolivet       if (!auglag->C) PetscCall(VecDuplicate(auglag->Y, &auglag->C));
1999371c9d4SSatish Balay     } else {
200ad540459SPierre Jolivet       if (!auglag->C) auglag->C = auglag->Ci;
201ad540459SPierre Jolivet       if (!auglag->Y) auglag->Y = auglag->Yi;
202661095bbSAlp Dener     }
203661095bbSAlp Dener   } else {
204ad540459SPierre Jolivet     if (!auglag->P) auglag->P = auglag->Px;
205ad540459SPierre Jolivet     if (!auglag->G) auglag->G = auglag->LgradX;
206ad540459SPierre Jolivet     if (!auglag->C) auglag->C = auglag->Ce;
207ad540459SPierre Jolivet     if (!auglag->Y) auglag->Y = auglag->Ye;
208661095bbSAlp Dener   }
209661095bbSAlp Dener   /* initialize parameters */
210661095bbSAlp Dener   if (auglag->type == TAO_ALMM_PHR) {
211661095bbSAlp Dener     auglag->mu_fac = 10.0;
212661095bbSAlp Dener     auglag->yi_min = 0.0;
213661095bbSAlp Dener     auglag->ytol0  = 0.5;
214661095bbSAlp Dener     auglag->gtol0  = tao->gatol;
215661095bbSAlp Dener     if (tao->gatol_changed && tao->catol_changed) {
2169566063dSJacob Faibussowitsch       PetscCall(PetscInfo(tao, "TAOALMM with PHR: different gradient and constraint tolerances are not supported, setting catol = gatol\n"));
217661095bbSAlp Dener       tao->catol = tao->gatol;
218661095bbSAlp Dener     }
219661095bbSAlp Dener   }
220661095bbSAlp Dener   /* set the Lagrangian formulation type for the subsolver */
221661095bbSAlp Dener   switch (auglag->type) {
2229371c9d4SSatish Balay   case TAO_ALMM_CLASSIC: auglag->sub_obj = TaoALMMComputeAugLagAndGradient_Private; break;
2239371c9d4SSatish Balay   case TAO_ALMM_PHR: auglag->sub_obj = TaoALMMComputePHRLagAndGradient_Private; break;
2249371c9d4SSatish Balay   default: break;
225661095bbSAlp Dener   }
226661095bbSAlp Dener   /* set up the subsolver */
2279566063dSJacob Faibussowitsch   PetscCall(TaoSetSolution(auglag->subsolver, auglag->P));
2289566063dSJacob Faibussowitsch   PetscCall(TaoSetObjective(auglag->subsolver, TaoALMMSubsolverObjective_Private, (void *)auglag));
2299566063dSJacob Faibussowitsch   PetscCall(TaoSetObjectiveAndGradient(auglag->subsolver, NULL, TaoALMMSubsolverObjectiveAndGradient_Private, (void *)auglag));
230661095bbSAlp Dener   if (tao->bounded) {
231661095bbSAlp Dener     /* make sure that the subsolver is a bound-constrained method */
2329566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)auglag->subsolver, TAOCG, &is_cg));
2339566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)auglag->subsolver, TAOLMVM, &is_lmvm));
234661095bbSAlp Dener     if (is_cg) {
2359566063dSJacob Faibussowitsch       PetscCall(TaoSetType(auglag->subsolver, TAOBNCG));
2369566063dSJacob Faibussowitsch       PetscCall(PetscInfo(tao, "TAOCG detected for bound-constrained problem, switching to TAOBNCG instead."));
237661095bbSAlp Dener     }
238661095bbSAlp Dener     if (is_lmvm) {
2399566063dSJacob Faibussowitsch       PetscCall(TaoSetType(auglag->subsolver, TAOBQNLS));
2409566063dSJacob Faibussowitsch       PetscCall(PetscInfo(tao, "TAOLMVM detected for bound-constrained problem, switching to TAOBQNLS instead."));
241661095bbSAlp Dener     }
242661095bbSAlp Dener     /* create lower and upper bound clone vectors for subsolver */
24348a46eb9SPierre Jolivet     if (!auglag->PL) PetscCall(VecDuplicate(auglag->P, &auglag->PL));
24448a46eb9SPierre Jolivet     if (!auglag->PU) PetscCall(VecDuplicate(auglag->P, &auglag->PU));
245661095bbSAlp Dener     if (tao->ineq_constrained) {
246661095bbSAlp Dener       /* create lower and upper bounds for slack, set lower to 0 */
2479566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(auglag->Ci, &SL));
2489566063dSJacob Faibussowitsch       PetscCall(VecSet(SL, 0.0));
2499566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(auglag->Ci, &SU));
2509566063dSJacob Faibussowitsch       PetscCall(VecSet(SU, PETSC_INFINITY));
251661095bbSAlp Dener       /* combine opt var bounds with slack bounds */
2529566063dSJacob Faibussowitsch       PetscCall(TaoALMMCombinePrimal_Private(tao, tao->XL, SL, auglag->PL));
2539566063dSJacob Faibussowitsch       PetscCall(TaoALMMCombinePrimal_Private(tao, tao->XU, SU, auglag->PU));
254661095bbSAlp Dener       /* destroy work vectors */
2559566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&SL));
2569566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&SU));
257661095bbSAlp Dener     } else {
258661095bbSAlp Dener       /* no inequality constraints, just copy bounds into the subsolver */
2599566063dSJacob Faibussowitsch       PetscCall(VecCopy(tao->XL, auglag->PL));
2609566063dSJacob Faibussowitsch       PetscCall(VecCopy(tao->XU, auglag->PU));
261661095bbSAlp Dener     }
2629566063dSJacob Faibussowitsch     PetscCall(TaoSetVariableBounds(auglag->subsolver, auglag->PL, auglag->PU));
263661095bbSAlp Dener   }
2649566063dSJacob Faibussowitsch   PetscCall(TaoSetUp(auglag->subsolver));
265661095bbSAlp Dener   auglag->G = auglag->subsolver->gradient;
266661095bbSAlp Dener 
267661095bbSAlp Dener   PetscFunctionReturn(0);
268661095bbSAlp Dener }
269661095bbSAlp Dener 
2709371c9d4SSatish Balay static PetscErrorCode TaoDestroy_ALMM(Tao tao) {
271661095bbSAlp Dener   TAO_ALMM *auglag = (TAO_ALMM *)tao->data;
272661095bbSAlp Dener 
273661095bbSAlp Dener   PetscFunctionBegin;
2749566063dSJacob Faibussowitsch   PetscCall(TaoDestroy(&auglag->subsolver));
275661095bbSAlp Dener   if (tao->setupcalled) {
2769566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&auglag->Xwork)); /* opt work */
277661095bbSAlp Dener     if (tao->eq_constrained) {
2789566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&auglag->Ye));     /* equality multipliers */
2799566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&auglag->Cework)); /* equality work vector */
280661095bbSAlp Dener     }
281661095bbSAlp Dener     if (tao->ineq_constrained) {
2829566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&auglag->Ps));     /* slack vars */
283661095bbSAlp Dener       auglag->Parr[0] = NULL;                 /* clear pointer to tao->solution, will be destroyed by TaoDestroy() shell */
2849566063dSJacob Faibussowitsch       PetscCall(PetscFree(auglag->Parr));     /* array of primal vectors */
2859566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&auglag->LgradS)); /* slack grad */
2869566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&auglag->Cizero)); /* zero vector for pointwise max */
2879566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&auglag->Yi));     /* inequality multipliers */
2889566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&auglag->Ciwork)); /* inequality work vector */
2899566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&auglag->P));      /* combo primal */
2909566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&auglag->Pis[0]));  /* index set for X inside P */
2919566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&auglag->Pis[1]));  /* index set for S inside P */
2929566063dSJacob Faibussowitsch       PetscCall(PetscFree(auglag->Pis));      /* array of P index sets */
2939566063dSJacob Faibussowitsch       PetscCall(VecScatterDestroy(&auglag->Pscatter[0]));
2949566063dSJacob Faibussowitsch       PetscCall(VecScatterDestroy(&auglag->Pscatter[1]));
2959566063dSJacob Faibussowitsch       PetscCall(PetscFree(auglag->Pscatter));
296661095bbSAlp Dener       if (tao->eq_constrained) {
2979566063dSJacob Faibussowitsch         PetscCall(VecDestroy(&auglag->Y));     /* combo multipliers */
2989566063dSJacob Faibussowitsch         PetscCall(PetscFree(auglag->Yarr));    /* array of dual vectors */
2999566063dSJacob Faibussowitsch         PetscCall(VecDestroy(&auglag->C));     /* combo constraints */
3009566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&auglag->Yis[0])); /* index set for Ye inside Y */
3019566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&auglag->Yis[1])); /* index set for Yi inside Y */
3029566063dSJacob Faibussowitsch         PetscCall(PetscFree(auglag->Yis));
3039566063dSJacob Faibussowitsch         PetscCall(VecScatterDestroy(&auglag->Yscatter[0]));
3049566063dSJacob Faibussowitsch         PetscCall(VecScatterDestroy(&auglag->Yscatter[1]));
3059566063dSJacob Faibussowitsch         PetscCall(PetscFree(auglag->Yscatter));
306661095bbSAlp Dener       }
307661095bbSAlp Dener     }
308661095bbSAlp Dener     if (tao->bounded) {
3099566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&auglag->PL)); /* lower bounds for subsolver */
3109566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&auglag->PU)); /* upper bounds for subsolver */
311661095bbSAlp Dener     }
312661095bbSAlp Dener   }
3132e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetType_C", NULL));
3142e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMSetType_C", NULL));
3152e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetSubsolver_C", NULL));
3162e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMSetSubsolver_C", NULL));
3172e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetMultipliers_C", NULL));
3182e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMSetMultipliers_C", NULL));
3192e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetPrimalIS_C", NULL));
3202e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetDualIS_C", NULL));
3219566063dSJacob Faibussowitsch   PetscCall(PetscFree(tao->data));
322661095bbSAlp Dener   PetscFunctionReturn(0);
323661095bbSAlp Dener }
324661095bbSAlp Dener 
3259371c9d4SSatish Balay static PetscErrorCode TaoSetFromOptions_ALMM(Tao tao, PetscOptionItems *PetscOptionsObject) {
326661095bbSAlp Dener   TAO_ALMM *auglag = (TAO_ALMM *)tao->data;
327661095bbSAlp Dener   PetscInt  i;
328661095bbSAlp Dener 
329661095bbSAlp Dener   PetscFunctionBegin;
33071075aafSSatish Balay   PetscOptionsHeadBegin(PetscOptionsObject, "Augmented Lagrangian multiplier method solves problems with general constraints by converting them into a sequence of unconstrained problems.");
3319566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_almm_mu_init", "initial penalty parameter", "", auglag->mu0, &auglag->mu0, NULL));
3329566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_almm_mu_factor", "increase factor for the penalty parameter", "", auglag->mu_fac, &auglag->mu_fac, NULL));
3339566063dSJacob Faibussowitsch   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));
3349566063dSJacob Faibussowitsch   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));
3359566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_almm_mu_max", "maximum safeguard for penalty parameter updates", "", auglag->mu_max, &auglag->mu_max, NULL));
3369566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_almm_ye_min", "minimum safeguard for equality multiplier updates", "", auglag->ye_min, &auglag->ye_min, NULL));
3379566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_almm_ye_max", "maximum safeguard for equality multipliers updates", "", auglag->ye_max, &auglag->ye_max, NULL));
3389566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_almm_yi_min", "minimum safeguard for inequality multipliers updates", "", auglag->yi_min, &auglag->yi_min, NULL));
3399566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_almm_yi_max", "maximum safeguard for inequality multipliers updates", "", auglag->yi_max, &auglag->yi_max, NULL));
3409566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-tao_almm_type", "augmented Lagrangian formulation type for the subproblem", "TaoALMMType", TaoALMMTypes, (PetscEnum)auglag->type, (PetscEnum *)&auglag->type, NULL));
341d0609cedSBarry Smith   PetscOptionsHeadEnd();
3429566063dSJacob Faibussowitsch   PetscCall(TaoSetOptionsPrefix(auglag->subsolver, ((PetscObject)tao)->prefix));
3439566063dSJacob Faibussowitsch   PetscCall(TaoAppendOptionsPrefix(auglag->subsolver, "tao_almm_subsolver_"));
3449566063dSJacob Faibussowitsch   PetscCall(TaoSetFromOptions(auglag->subsolver));
345661095bbSAlp Dener   for (i = 0; i < tao->numbermonitors; i++) {
3469566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)tao->monitorcontext[i]));
3479566063dSJacob Faibussowitsch     PetscCall(TaoSetMonitor(auglag->subsolver, tao->monitor[i], tao->monitorcontext[i], tao->monitordestroy[i]));
348ad540459SPierre Jolivet     if (tao->monitor[i] == TaoMonitorDefault || tao->monitor[i] == TaoDefaultCMonitor || tao->monitor[i] == TaoDefaultGMonitor || tao->monitor[i] == TaoDefaultSMonitor) auglag->info = PETSC_TRUE;
349661095bbSAlp Dener   }
350661095bbSAlp Dener   PetscFunctionReturn(0);
351661095bbSAlp Dener }
352661095bbSAlp Dener 
353661095bbSAlp Dener /* -------------------------------------------------------- */
354661095bbSAlp Dener 
355661095bbSAlp Dener /*MC
356661095bbSAlp Dener   TaoALMM - Augmented Lagrangian multiplier method for solving nonlinear optimization problems with general constraints.
357661095bbSAlp Dener 
358661095bbSAlp Dener   Options Database Keys:
359661095bbSAlp Dener + -tao_almm_mu_init <real>       - initial penalty parameter (default: 10.)
360661095bbSAlp Dener . -tao_almm_mu_factor <real>     - increase factor for the penalty parameter (default: 100.)
361661095bbSAlp Dener . -tao_almm_mu_max <real>        - maximum safeguard for penalty parameter updates (default: 1.e20)
362661095bbSAlp Dener . -tao_almm_mu_power_good <real> - exponential for penalty parameter when multiplier update is accepted (default: 0.9)
363661095bbSAlp Dener . -tao_almm_mu_power_bad <real>  - exponential for penalty parameter when multiplier update is rejected (default: 0.1)
364661095bbSAlp Dener . -tao_almm_ye_min <real>        - minimum safeguard for equality multiplier updates (default: -1.e20)
365661095bbSAlp Dener . -tao_almm_ye_max <real>        - maximum safeguard for equality multiplier updates (default: 1.e20)
366661095bbSAlp Dener . -tao_almm_yi_min <real>        - minimum safeguard for inequality multiplier updates (default: -1.e20)
367661095bbSAlp Dener . -tao_almm_yi_max <real>        - maximum safeguard for inequality multiplier updates (default: 1.e20)
368661095bbSAlp Dener - -tao_almm_type <classic,phr>   - change formulation of the augmented Lagrangian merit function for the subproblem (default: classic)
369661095bbSAlp Dener 
370661095bbSAlp Dener   Level: beginner
371661095bbSAlp Dener 
372661095bbSAlp Dener   Notes:
373661095bbSAlp Dener   This method converts a constrained problem into a sequence of unconstrained problems via the augmented
374661095bbSAlp Dener   Lagrangian merit function. Bound constraints are pushed down to the subproblem without any modifications.
375661095bbSAlp Dener 
376661095bbSAlp Dener   Two formulations are offered for the subproblem: canonical Hestenes-Powell augmented Lagrangian with slack
377661095bbSAlp Dener   variables for inequality constraints, and a slack-less Powell-Hestenes-Rockafellar (PHR) formulation utilizing a
378661095bbSAlp Dener   pointwise max() penalty on inequality constraints. The canonical augmented Lagrangian formulation typically
379661095bbSAlp Dener   converges faster for most problems. However, PHR may be desirable for problems featuring a large number
380661095bbSAlp Dener   of inequality constraints because it avoids inflating the size of the subproblem with slack variables.
381661095bbSAlp Dener 
382661095bbSAlp Dener   The subproblem is solved using a nested first-order TAO solver. The user can retrieve a pointer to
383661095bbSAlp Dener   the subsolver via TaoALMMGetSubsolver() or pass command line arguments to it using the
384661095bbSAlp Dener   "-tao_almm_subsolver_" prefix. Currently, TaoALMM does not support second-order methods for the
385661095bbSAlp Dener   subproblem. It is also highly recommended that the subsolver chosen by the user utilize a trust-region
386661095bbSAlp Dener   strategy for globalization (default: TAOBQNKTR) especially if the outer problem features bound constraints.
387661095bbSAlp Dener 
388661095bbSAlp Dener .vb
389661095bbSAlp Dener   while unconverged
390661095bbSAlp Dener     solve argmin_x L(x) s.t. l <= x <= u
391661095bbSAlp Dener     if ||c|| <= y_tol
392661095bbSAlp Dener       if ||c|| <= c_tol && ||Lgrad|| <= g_tol:
393661095bbSAlp Dener         problem converged, return solution
394661095bbSAlp Dener       else
395661095bbSAlp Dener         constraints sufficiently improved
396661095bbSAlp Dener         update multipliers and tighten tolerances
397661095bbSAlp Dener       endif
398661095bbSAlp Dener     else
399661095bbSAlp Dener       constraints did not improve
400661095bbSAlp Dener       update penalty and loosen tolerances
401661095bbSAlp Dener     endif
402661095bbSAlp Dener   endwhile
403661095bbSAlp Dener .ve
404661095bbSAlp Dener 
405db781477SPatrick Sanan .seealso: `TaoALMMGetType()`, `TaoALMMSetType()`, `TaoALMMSetSubsolver()`, `TaoALMMGetSubsolver()`,
406db781477SPatrick Sanan           `TaoALMMGetMultipliers()`, `TaoALMMSetMultipliers()`, `TaoALMMGetPrimalIS()`, `TaoALMMGetDualIS()`
407661095bbSAlp Dener M*/
4089371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode TaoCreate_ALMM(Tao tao) {
409661095bbSAlp Dener   TAO_ALMM *auglag;
410661095bbSAlp Dener 
411661095bbSAlp Dener   PetscFunctionBegin;
412*4dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&auglag));
413661095bbSAlp Dener 
414661095bbSAlp Dener   tao->ops->destroy        = TaoDestroy_ALMM;
415661095bbSAlp Dener   tao->ops->setup          = TaoSetUp_ALMM;
416661095bbSAlp Dener   tao->ops->setfromoptions = TaoSetFromOptions_ALMM;
417661095bbSAlp Dener   tao->ops->view           = TaoView_ALMM;
418661095bbSAlp Dener   tao->ops->solve          = TaoSolve_ALMM;
419661095bbSAlp Dener 
420661095bbSAlp Dener   tao->gatol = 1.e-5;
421661095bbSAlp Dener   tao->grtol = 0.0;
422661095bbSAlp Dener   tao->gttol = 0.0;
423661095bbSAlp Dener   tao->catol = 1.e-5;
424661095bbSAlp Dener   tao->crtol = 0.0;
425661095bbSAlp Dener 
426661095bbSAlp Dener   tao->data           = (void *)auglag;
427661095bbSAlp Dener   auglag->parent      = tao;
428661095bbSAlp Dener   auglag->mu0         = 10.0;
429661095bbSAlp Dener   auglag->mu          = auglag->mu0;
430661095bbSAlp Dener   auglag->mu_fac      = 10.0;
431661095bbSAlp Dener   auglag->mu_max      = PETSC_INFINITY;
432661095bbSAlp Dener   auglag->mu_pow_good = 0.9;
433661095bbSAlp Dener   auglag->mu_pow_bad  = 0.1;
434661095bbSAlp Dener   auglag->ye_min      = PETSC_NINFINITY;
435661095bbSAlp Dener   auglag->ye_max      = PETSC_INFINITY;
436661095bbSAlp Dener   auglag->yi_min      = PETSC_NINFINITY;
437661095bbSAlp Dener   auglag->yi_max      = PETSC_INFINITY;
438661095bbSAlp Dener   auglag->ytol0       = 0.1 / PetscPowReal(auglag->mu0, auglag->mu_pow_bad);
439661095bbSAlp Dener   auglag->ytol        = auglag->ytol0;
440661095bbSAlp Dener   auglag->gtol0       = 1.0 / auglag->mu0;
441661095bbSAlp Dener   auglag->gtol        = auglag->gtol0;
442661095bbSAlp Dener 
443661095bbSAlp Dener   auglag->sub_obj = TaoALMMComputeAugLagAndGradient_Private;
444661095bbSAlp Dener   auglag->type    = TAO_ALMM_CLASSIC;
445661095bbSAlp Dener   auglag->info    = PETSC_FALSE;
446661095bbSAlp Dener 
4479566063dSJacob Faibussowitsch   PetscCall(TaoCreate(PetscObjectComm((PetscObject)tao), &auglag->subsolver));
4489566063dSJacob Faibussowitsch   PetscCall(TaoSetType(auglag->subsolver, TAOBQNKTR));
4499566063dSJacob Faibussowitsch   PetscCall(TaoSetTolerances(auglag->subsolver, auglag->gtol, 0.0, 0.0));
4509566063dSJacob Faibussowitsch   PetscCall(TaoSetMaximumIterations(auglag->subsolver, 1000));
4519566063dSJacob Faibussowitsch   PetscCall(TaoSetMaximumFunctionEvaluations(auglag->subsolver, 10000));
4529566063dSJacob Faibussowitsch   PetscCall(TaoSetFunctionLowerBound(auglag->subsolver, PETSC_NINFINITY));
4539566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)auglag->subsolver, (PetscObject)tao, 1));
454661095bbSAlp Dener 
4559566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetType_C", TaoALMMGetType_Private));
4569566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMSetType_C", TaoALMMSetType_Private));
4579566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetSubsolver_C", TaoALMMGetSubsolver_Private));
4589566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMSetSubsolver_C", TaoALMMSetSubsolver_Private));
4599566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetMultipliers_C", TaoALMMGetMultipliers_Private));
4609566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMSetMultipliers_C", TaoALMMSetMultipliers_Private));
4619566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetPrimalIS_C", TaoALMMGetPrimalIS_Private));
4629566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetDualIS_C", TaoALMMGetDualIS_Private));
463661095bbSAlp Dener   PetscFunctionReturn(0);
464661095bbSAlp Dener }
465661095bbSAlp Dener 
4669371c9d4SSatish Balay static PetscErrorCode TaoALMMCombinePrimal_Private(Tao tao, Vec X, Vec S, Vec P) {
467661095bbSAlp Dener   TAO_ALMM *auglag = (TAO_ALMM *)tao->data;
468661095bbSAlp Dener 
469661095bbSAlp Dener   PetscFunctionBegin;
470661095bbSAlp Dener   if (tao->ineq_constrained) {
4719566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(auglag->Pscatter[0], X, P, INSERT_VALUES, SCATTER_REVERSE));
4729566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(auglag->Pscatter[0], X, P, INSERT_VALUES, SCATTER_REVERSE));
4739566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(auglag->Pscatter[1], S, P, INSERT_VALUES, SCATTER_REVERSE));
4749566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(auglag->Pscatter[1], S, P, INSERT_VALUES, SCATTER_REVERSE));
475661095bbSAlp Dener   } else {
4769566063dSJacob Faibussowitsch     PetscCall(VecCopy(X, P));
477661095bbSAlp Dener   }
478661095bbSAlp Dener   PetscFunctionReturn(0);
479661095bbSAlp Dener }
480661095bbSAlp Dener 
4819371c9d4SSatish Balay static PetscErrorCode TaoALMMCombineDual_Private(Tao tao, Vec EQ, Vec IN, Vec Y) {
482661095bbSAlp Dener   TAO_ALMM *auglag = (TAO_ALMM *)tao->data;
483661095bbSAlp Dener 
484661095bbSAlp Dener   PetscFunctionBegin;
485661095bbSAlp Dener   if (tao->eq_constrained) {
486661095bbSAlp Dener     if (tao->ineq_constrained) {
4879566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(auglag->Yscatter[0], EQ, Y, INSERT_VALUES, SCATTER_REVERSE));
4889566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(auglag->Yscatter[0], EQ, Y, INSERT_VALUES, SCATTER_REVERSE));
4899566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(auglag->Yscatter[1], IN, Y, INSERT_VALUES, SCATTER_REVERSE));
4909566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(auglag->Yscatter[1], IN, Y, INSERT_VALUES, SCATTER_REVERSE));
491661095bbSAlp Dener     } else {
4929566063dSJacob Faibussowitsch       PetscCall(VecCopy(EQ, Y));
493661095bbSAlp Dener     }
494661095bbSAlp Dener   } else {
4959566063dSJacob Faibussowitsch     PetscCall(VecCopy(IN, Y));
496661095bbSAlp Dener   }
497661095bbSAlp Dener   PetscFunctionReturn(0);
498661095bbSAlp Dener }
499661095bbSAlp Dener 
5009371c9d4SSatish Balay static PetscErrorCode TaoALMMSplitPrimal_Private(Tao tao, Vec P, Vec X, Vec S) {
501661095bbSAlp Dener   TAO_ALMM *auglag = (TAO_ALMM *)tao->data;
502661095bbSAlp Dener 
503661095bbSAlp Dener   PetscFunctionBegin;
504661095bbSAlp Dener   if (tao->ineq_constrained) {
5059566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(auglag->Pscatter[0], P, X, INSERT_VALUES, SCATTER_FORWARD));
5069566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(auglag->Pscatter[0], P, X, INSERT_VALUES, SCATTER_FORWARD));
5079566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(auglag->Pscatter[1], P, S, INSERT_VALUES, SCATTER_FORWARD));
5089566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(auglag->Pscatter[1], P, S, INSERT_VALUES, SCATTER_FORWARD));
509661095bbSAlp Dener   } else {
5109566063dSJacob Faibussowitsch     PetscCall(VecCopy(P, X));
511661095bbSAlp Dener   }
512661095bbSAlp Dener   PetscFunctionReturn(0);
513661095bbSAlp Dener }
514661095bbSAlp Dener 
515661095bbSAlp Dener /* this assumes that the latest constraints are stored in Ce and Ci, and also combined in C */
5169371c9d4SSatish Balay static PetscErrorCode TaoALMMComputeOptimalityNorms_Private(Tao tao) {
517661095bbSAlp Dener   TAO_ALMM *auglag = (TAO_ALMM *)tao->data;
518661095bbSAlp Dener 
519661095bbSAlp Dener   PetscFunctionBegin;
520661095bbSAlp Dener   /* if bounded, project the gradient */
5211baa6e33SBarry Smith   if (tao->bounded) PetscCall(VecBoundGradientProjection(auglag->LgradX, auglag->Px, tao->XL, tao->XU, auglag->LgradX));
522661095bbSAlp Dener   if (auglag->type == TAO_ALMM_PHR) {
5239566063dSJacob Faibussowitsch     PetscCall(VecNorm(auglag->LgradX, NORM_INFINITY, &auglag->gnorm));
524661095bbSAlp Dener     auglag->cenorm = 0.0;
52548a46eb9SPierre Jolivet     if (tao->eq_constrained) PetscCall(VecNorm(auglag->Ce, NORM_INFINITY, &auglag->cenorm));
526661095bbSAlp Dener     auglag->cinorm = 0.0;
527661095bbSAlp Dener     if (tao->ineq_constrained) {
5289566063dSJacob Faibussowitsch       PetscCall(VecCopy(auglag->Yi, auglag->Ciwork));
5299566063dSJacob Faibussowitsch       PetscCall(VecScale(auglag->Ciwork, -1.0 / auglag->mu));
5309566063dSJacob Faibussowitsch       PetscCall(VecPointwiseMax(auglag->Ciwork, auglag->Ci, auglag->Ciwork));
5319566063dSJacob Faibussowitsch       PetscCall(VecNorm(auglag->Ciwork, NORM_INFINITY, &auglag->cinorm));
532661095bbSAlp Dener     }
533661095bbSAlp Dener     auglag->cnorm_old = auglag->cnorm;
534661095bbSAlp Dener     auglag->cnorm     = PetscMax(auglag->cenorm, auglag->cinorm);
535661095bbSAlp Dener     auglag->ytol      = auglag->ytol0 * auglag->cnorm_old;
536661095bbSAlp Dener   } else {
5379566063dSJacob Faibussowitsch     PetscCall(VecNorm(auglag->LgradX, NORM_2, &auglag->gnorm));
5389566063dSJacob Faibussowitsch     PetscCall(VecNorm(auglag->C, NORM_2, &auglag->cnorm));
539661095bbSAlp Dener   }
540661095bbSAlp Dener   PetscFunctionReturn(0);
541661095bbSAlp Dener }
542661095bbSAlp Dener 
5439371c9d4SSatish Balay static PetscErrorCode TaoALMMEvaluateIterate_Private(Tao tao, Vec P) {
544661095bbSAlp Dener   TAO_ALMM *auglag = (TAO_ALMM *)tao->data;
545661095bbSAlp Dener 
546661095bbSAlp Dener   PetscFunctionBegin;
547661095bbSAlp Dener   /* split solution into primal and slack components */
5489566063dSJacob Faibussowitsch   PetscCall(TaoALMMSplitPrimal_Private(tao, auglag->P, auglag->Px, auglag->Ps));
549661095bbSAlp Dener 
550661095bbSAlp Dener   /* compute f, df/dx and the constraints */
5519566063dSJacob Faibussowitsch   PetscCall(TaoComputeObjectiveAndGradient(tao, auglag->Px, &auglag->fval, auglag->LgradX));
552661095bbSAlp Dener   if (tao->eq_constrained) {
5539566063dSJacob Faibussowitsch     PetscCall(TaoComputeEqualityConstraints(tao, auglag->Px, auglag->Ce));
5549566063dSJacob Faibussowitsch     PetscCall(TaoComputeJacobianEquality(tao, auglag->Px, auglag->Ae, auglag->Ae));
555661095bbSAlp Dener   }
556661095bbSAlp Dener   if (tao->ineq_constrained) {
5579566063dSJacob Faibussowitsch     PetscCall(TaoComputeInequalityConstraints(tao, auglag->Px, auglag->Ci));
5589566063dSJacob Faibussowitsch     PetscCall(TaoComputeJacobianInequality(tao, auglag->Px, auglag->Ai, auglag->Ai));
559661095bbSAlp Dener     switch (auglag->type) {
5607721a36fSStefano Zampini     case TAO_ALMM_CLASSIC:
561661095bbSAlp Dener       /* classic formulation converts inequality to equality constraints via slack variables */
5629566063dSJacob Faibussowitsch       PetscCall(VecAXPY(auglag->Ci, -1.0, auglag->Ps));
563661095bbSAlp Dener       break;
5647721a36fSStefano Zampini     case TAO_ALMM_PHR:
565661095bbSAlp Dener       /* PHR is based on Ci <= 0 while TAO defines Ci >= 0 so we hit it with a negative sign */
5669566063dSJacob Faibussowitsch       PetscCall(VecScale(auglag->Ci, -1.0));
5679566063dSJacob Faibussowitsch       PetscCall(MatScale(auglag->Ai, -1.0));
568661095bbSAlp Dener       break;
5699371c9d4SSatish Balay     default: break;
570661095bbSAlp Dener     }
571661095bbSAlp Dener   }
572661095bbSAlp Dener   /* combine constraints into one vector */
5739566063dSJacob Faibussowitsch   PetscCall(TaoALMMCombineDual_Private(tao, auglag->Ce, auglag->Ci, auglag->C));
574661095bbSAlp Dener   PetscFunctionReturn(0);
575661095bbSAlp Dener }
576661095bbSAlp Dener 
577661095bbSAlp Dener /*
578661095bbSAlp Dener Lphr = f + 0.5*mu*[ (Ce + Ye/mu)^T (Ce + Ye/mu) + pmin(0, Ci + Yi/mu)^T pmin(0, Ci + Yi/mu)]
579661095bbSAlp Dener 
580661095bbSAlp Dener dLphr/dX = dF/dX + mu*[ (Ce + Ye/mu)^T Ae + pmin(0, Ci + Yi/mu)^T Ai]
581661095bbSAlp Dener 
582661095bbSAlp Dener dLphr/dS = 0
583661095bbSAlp Dener */
5849371c9d4SSatish Balay static PetscErrorCode TaoALMMComputePHRLagAndGradient_Private(Tao tao) {
585661095bbSAlp Dener   TAO_ALMM *auglag  = (TAO_ALMM *)tao->data;
586661095bbSAlp Dener   PetscReal eq_norm = 0.0, ineq_norm = 0.0;
587661095bbSAlp Dener 
588661095bbSAlp Dener   PetscFunctionBegin;
5899566063dSJacob Faibussowitsch   PetscCall(TaoALMMEvaluateIterate_Private(tao, auglag->P));
590661095bbSAlp Dener   if (tao->eq_constrained) {
591661095bbSAlp Dener     /* Ce_work = mu*(Ce + Ye/mu) */
5929566063dSJacob Faibussowitsch     PetscCall(VecWAXPY(auglag->Cework, 1.0 / auglag->mu, auglag->Ye, auglag->Ce));
5939566063dSJacob Faibussowitsch     PetscCall(VecDot(auglag->Cework, auglag->Cework, &eq_norm)); /* contribution to scalar Lagrangian */
5949566063dSJacob Faibussowitsch     PetscCall(VecScale(auglag->Cework, auglag->mu));
595661095bbSAlp Dener     /* dL/dX += mu*(Ce + Ye/mu)^T Ae */
5969566063dSJacob Faibussowitsch     PetscCall(MatMultTransposeAdd(auglag->Ae, auglag->Cework, auglag->LgradX, auglag->LgradX));
597661095bbSAlp Dener   }
598661095bbSAlp Dener   if (tao->ineq_constrained) {
599661095bbSAlp Dener     /* Ci_work = mu * pmax(0, Ci + Yi/mu) where pmax() is pointwise max() */
6009566063dSJacob Faibussowitsch     PetscCall(VecWAXPY(auglag->Ciwork, 1.0 / auglag->mu, auglag->Yi, auglag->Ci));
6019566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMax(auglag->Ciwork, auglag->Cizero, auglag->Ciwork));
6029566063dSJacob Faibussowitsch     PetscCall(VecDot(auglag->Ciwork, auglag->Ciwork, &ineq_norm)); /* contribution to scalar Lagrangian */
603661095bbSAlp Dener     /* dL/dX += mu * pmax(0, Ci + Yi/mu)^T Ai */
6049566063dSJacob Faibussowitsch     PetscCall(VecScale(auglag->Ciwork, auglag->mu));
6059566063dSJacob Faibussowitsch     PetscCall(MatMultTransposeAdd(auglag->Ai, auglag->Ciwork, auglag->LgradX, auglag->LgradX));
606a5b23f4aSJose E. Roman     /* dL/dS = 0 because there are no slacks in PHR */
6079566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(auglag->LgradS));
608661095bbSAlp Dener   }
609661095bbSAlp Dener   /* combine gradient together */
6109566063dSJacob Faibussowitsch   PetscCall(TaoALMMCombinePrimal_Private(tao, auglag->LgradX, auglag->LgradS, auglag->G));
611661095bbSAlp Dener   /* compute L = f + 0.5 * mu * [(Ce + Ye/mu)^T (Ce + Ye/mu) + pmax(0, Ci + Yi/mu)^T pmax(0, Ci + Yi/mu)] */
6125f80ce2aSJacob Faibussowitsch   auglag->Lval = auglag->fval + 0.5 * auglag->mu * (eq_norm + ineq_norm);
613661095bbSAlp Dener   PetscFunctionReturn(0);
614661095bbSAlp Dener }
615661095bbSAlp Dener 
616661095bbSAlp Dener /*
617661095bbSAlp Dener Lc = F + Ye^TCe + Yi^T(Ci - S) + 0.5*mu*[Ce^TCe + (Ci - S)^T(Ci - S)]
618661095bbSAlp Dener 
619661095bbSAlp Dener dLc/dX = dF/dX + Ye^TAe + Yi^TAi + 0.5*mu*[Ce^TAe + (Ci - S)^TAi]
620661095bbSAlp Dener 
621661095bbSAlp Dener dLc/dS = -[Yi + mu*(Ci - S)]
622661095bbSAlp Dener */
6239371c9d4SSatish Balay static PetscErrorCode TaoALMMComputeAugLagAndGradient_Private(Tao tao) {
624661095bbSAlp Dener   TAO_ALMM *auglag = (TAO_ALMM *)tao->data;
625661095bbSAlp Dener   PetscReal yeTce = 0.0, yiTcims = 0.0, ceTce = 0.0, cimsTcims = 0.0;
626661095bbSAlp Dener 
627661095bbSAlp Dener   PetscFunctionBegin;
6289566063dSJacob Faibussowitsch   PetscCall(TaoALMMEvaluateIterate_Private(tao, auglag->P));
629661095bbSAlp Dener   if (tao->eq_constrained) {
630661095bbSAlp Dener     /* compute scalar contributions */
6319566063dSJacob Faibussowitsch     PetscCall(VecDot(auglag->Ye, auglag->Ce, &yeTce));
6329566063dSJacob Faibussowitsch     PetscCall(VecDot(auglag->Ce, auglag->Ce, &ceTce));
633661095bbSAlp Dener     /* dL/dX += ye^T Ae */
6349566063dSJacob Faibussowitsch     PetscCall(MatMultTransposeAdd(auglag->Ae, auglag->Ye, auglag->LgradX, auglag->LgradX));
635661095bbSAlp Dener     /* dL/dX += mu * ce^T Ae */
6369566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(auglag->Ae, auglag->Ce, auglag->Xwork));
6379566063dSJacob Faibussowitsch     PetscCall(VecAXPY(auglag->LgradX, auglag->mu, auglag->Xwork));
638661095bbSAlp Dener   }
639661095bbSAlp Dener   if (tao->ineq_constrained) {
640661095bbSAlp Dener     /* compute scalar contributions */
6419566063dSJacob Faibussowitsch     PetscCall(VecDot(auglag->Yi, auglag->Ci, &yiTcims));
6429566063dSJacob Faibussowitsch     PetscCall(VecDot(auglag->Ci, auglag->Ci, &cimsTcims));
643661095bbSAlp Dener     /* dL/dX += yi^T Ai */
6449566063dSJacob Faibussowitsch     PetscCall(MatMultTransposeAdd(auglag->Ai, auglag->Yi, auglag->LgradX, auglag->LgradX));
645661095bbSAlp Dener     /* dL/dX += mu * (ci - s)^T Ai */
6469566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(auglag->Ai, auglag->Ci, auglag->Xwork));
6479566063dSJacob Faibussowitsch     PetscCall(VecAXPY(auglag->LgradX, auglag->mu, auglag->Xwork));
648661095bbSAlp Dener     /* dL/dS = -[yi + mu*(ci - s)] */
6499566063dSJacob Faibussowitsch     PetscCall(VecWAXPY(auglag->LgradS, auglag->mu, auglag->Ci, auglag->Yi));
6509566063dSJacob Faibussowitsch     PetscCall(VecScale(auglag->LgradS, -1.0));
651661095bbSAlp Dener   }
652661095bbSAlp Dener   /* combine gradient together */
6539566063dSJacob Faibussowitsch   PetscCall(TaoALMMCombinePrimal_Private(tao, auglag->LgradX, auglag->LgradS, auglag->G));
654661095bbSAlp Dener   /* compute L = f + ye^T ce + yi^T (ci - s) + 0.5*mu*||ce||^2 + 0.5*mu*||ci - s||^2 */
655661095bbSAlp Dener   auglag->Lval = auglag->fval + yeTce + yiTcims + 0.5 * auglag->mu * (ceTce + cimsTcims);
656661095bbSAlp Dener   PetscFunctionReturn(0);
657661095bbSAlp Dener }
658661095bbSAlp Dener 
6599371c9d4SSatish Balay PetscErrorCode TaoALMMSubsolverObjective_Private(Tao tao, Vec P, PetscReal *Lval, void *ctx) {
6607721a36fSStefano Zampini   TAO_ALMM *auglag = (TAO_ALMM *)ctx;
6617721a36fSStefano Zampini 
6627721a36fSStefano Zampini   PetscFunctionBegin;
6639566063dSJacob Faibussowitsch   PetscCall(VecCopy(P, auglag->P));
6649566063dSJacob Faibussowitsch   PetscCall((*auglag->sub_obj)(auglag->parent));
6657721a36fSStefano Zampini   *Lval = auglag->Lval;
6667721a36fSStefano Zampini   PetscFunctionReturn(0);
6677721a36fSStefano Zampini }
6687721a36fSStefano Zampini 
6699371c9d4SSatish Balay PetscErrorCode TaoALMMSubsolverObjectiveAndGradient_Private(Tao tao, Vec P, PetscReal *Lval, Vec G, void *ctx) {
670661095bbSAlp Dener   TAO_ALMM *auglag = (TAO_ALMM *)ctx;
671661095bbSAlp Dener 
672661095bbSAlp Dener   PetscFunctionBegin;
6739566063dSJacob Faibussowitsch   PetscCall(VecCopy(P, auglag->P));
6749566063dSJacob Faibussowitsch   PetscCall((*auglag->sub_obj)(auglag->parent));
6759566063dSJacob Faibussowitsch   PetscCall(VecCopy(auglag->G, G));
676661095bbSAlp Dener   *Lval = auglag->Lval;
677661095bbSAlp Dener   PetscFunctionReturn(0);
678661095bbSAlp Dener }
679