xref: /petsc/src/tao/constrained/impls/almm/almm.c (revision c8b0c2b58db46088bcb0f3a2497e96cb9efa794b)
1954e39ddSJose E. Roman #include <../src/tao/constrained/impls/almm/almm.h> /*I "petsctao.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 
13d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSolve_ALMM(Tao tao)
14d71ae5a4SJacob Faibussowitsch {
15661095bbSAlp Dener   TAO_ALMM          *auglag = (TAO_ALMM *)tao->data;
16661095bbSAlp Dener   TaoConvergedReason reason;
17661095bbSAlp Dener   PetscReal          updated;
18661095bbSAlp Dener 
19661095bbSAlp Dener   PetscFunctionBegin;
20661095bbSAlp Dener   /* reset initial multiplier/slack guess */
21661095bbSAlp Dener   if (!tao->recycle) {
22661095bbSAlp Dener     if (tao->ineq_constrained) {
239566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(auglag->Ps));
249566063dSJacob Faibussowitsch       PetscCall(TaoALMMCombinePrimal_Private(tao, auglag->Px, auglag->Ps, auglag->P));
25*c8b0c2b5SHansol Suh       PetscCall(VecSet(auglag->Yi, 0.0));
26661095bbSAlp Dener     }
27*c8b0c2b5SHansol Suh     if (tao->eq_constrained) PetscCall(VecSet(auglag->Ye, 0.0));
28661095bbSAlp Dener   }
29661095bbSAlp Dener 
30661095bbSAlp Dener   /* compute initial nonlinear Lagrangian and its derivatives */
319566063dSJacob Faibussowitsch   PetscCall((*auglag->sub_obj)(tao));
329566063dSJacob Faibussowitsch   PetscCall(TaoALMMComputeOptimalityNorms_Private(tao));
33661095bbSAlp Dener   /* print initial step and check convergence */
349566063dSJacob Faibussowitsch   PetscCall(PetscInfo(tao, "Solving with %s formulation\n", TaoALMMTypes[auglag->type]));
359566063dSJacob Faibussowitsch   PetscCall(TaoLogConvergenceHistory(tao, auglag->Lval, auglag->gnorm, auglag->cnorm, tao->ksp_its));
369566063dSJacob Faibussowitsch   PetscCall(TaoMonitor(tao, tao->niter, auglag->fval, auglag->gnorm, auglag->cnorm, 0.0));
37dbbe0bcdSBarry Smith   PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
38661095bbSAlp Dener   /* set initial penalty factor and inner solver tolerance */
39661095bbSAlp Dener   switch (auglag->type) {
40d71ae5a4SJacob Faibussowitsch   case TAO_ALMM_CLASSIC:
41d71ae5a4SJacob Faibussowitsch     auglag->mu = auglag->mu0;
42d71ae5a4SJacob Faibussowitsch     break;
437721a36fSStefano Zampini   case TAO_ALMM_PHR:
44661095bbSAlp Dener     auglag->cenorm = 0.0;
4548a46eb9SPierre Jolivet     if (tao->eq_constrained) PetscCall(VecDot(auglag->Ce, auglag->Ce, &auglag->cenorm));
46661095bbSAlp Dener     auglag->cinorm = 0.0;
47661095bbSAlp Dener     if (tao->ineq_constrained) {
489566063dSJacob Faibussowitsch       PetscCall(VecCopy(auglag->Ci, auglag->Ciwork));
499566063dSJacob Faibussowitsch       PetscCall(VecScale(auglag->Ciwork, -1.0));
509566063dSJacob Faibussowitsch       PetscCall(VecPointwiseMax(auglag->Ciwork, auglag->Cizero, auglag->Ciwork));
519566063dSJacob Faibussowitsch       PetscCall(VecDot(auglag->Ciwork, auglag->Ciwork, &auglag->cinorm));
52661095bbSAlp Dener     }
53661095bbSAlp Dener     /* determine initial penalty factor based on the balance of constraint violation and objective function value */
5408fc2fa2SPaul T. Kühner     if (PetscAbsReal(auglag->cenorm + auglag->cinorm) == 0.0) auglag->mu = 10.0;
5508fc2fa2SPaul T. Kühner     else auglag->mu = PetscMax(1.e-6, PetscMin(10.0, 2.0 * PetscAbsReal(auglag->fval) / (auglag->cenorm + auglag->cinorm)));
56661095bbSAlp Dener     break;
57d71ae5a4SJacob Faibussowitsch   default:
58d71ae5a4SJacob Faibussowitsch     break;
59661095bbSAlp Dener   }
60661095bbSAlp Dener   auglag->gtol = auglag->gtol0;
6163a3b9bcSJacob Faibussowitsch   PetscCall(PetscInfo(tao, "Initial penalty: %.2g\n", (double)auglag->mu));
62661095bbSAlp Dener 
63661095bbSAlp Dener   /* start aug-lag outer loop */
64661095bbSAlp Dener   while (tao->reason == TAO_CONTINUE_ITERATING) {
65661095bbSAlp Dener     ++tao->niter;
66661095bbSAlp Dener     /* update subsolver tolerance */
6763a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(tao, "Subsolver tolerance: ||G|| <= %e\n", (double)auglag->gtol));
689566063dSJacob Faibussowitsch     PetscCall(TaoSetTolerances(auglag->subsolver, auglag->gtol, 0.0, 0.0));
69661095bbSAlp Dener     /* solve the bound-constrained or unconstrained subproblem */
709a09327dSAlp Dener     PetscCall(VecCopy(auglag->P, auglag->subsolver->solution));
719566063dSJacob Faibussowitsch     PetscCall(TaoSolve(auglag->subsolver));
729a09327dSAlp Dener     PetscCall(VecCopy(auglag->subsolver->solution, auglag->P));
739566063dSJacob Faibussowitsch     PetscCall(TaoGetConvergedReason(auglag->subsolver, &reason));
74661095bbSAlp Dener     tao->ksp_its += auglag->subsolver->ksp_its;
7548a46eb9SPierre Jolivet     if (reason != TAO_CONVERGED_GATOL) PetscCall(PetscInfo(tao, "Subsolver failed to converge, reason: %s\n", TaoConvergedReasons[reason]));
76661095bbSAlp Dener     /* evaluate solution and test convergence */
779566063dSJacob Faibussowitsch     PetscCall((*auglag->sub_obj)(tao));
789566063dSJacob Faibussowitsch     PetscCall(TaoALMMComputeOptimalityNorms_Private(tao));
79661095bbSAlp Dener     /* decide whether to update multipliers or not */
80661095bbSAlp Dener     updated = 0.0;
81661095bbSAlp Dener     if (auglag->cnorm <= auglag->ytol) {
8263a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(tao, "Multipliers updated: ||C|| <= %e\n", (double)auglag->ytol));
83661095bbSAlp Dener       /* constraints are good, update multipliers and convergence tolerances */
84661095bbSAlp Dener       if (tao->eq_constrained) {
859566063dSJacob Faibussowitsch         PetscCall(VecAXPY(auglag->Ye, auglag->mu, auglag->Ce));
869566063dSJacob Faibussowitsch         PetscCall(VecSet(auglag->Cework, auglag->ye_max));
879566063dSJacob Faibussowitsch         PetscCall(VecPointwiseMin(auglag->Ye, auglag->Cework, auglag->Ye));
889566063dSJacob Faibussowitsch         PetscCall(VecSet(auglag->Cework, auglag->ye_min));
899566063dSJacob Faibussowitsch         PetscCall(VecPointwiseMax(auglag->Ye, auglag->Cework, auglag->Ye));
90661095bbSAlp Dener       }
91661095bbSAlp Dener       if (tao->ineq_constrained) {
929566063dSJacob Faibussowitsch         PetscCall(VecAXPY(auglag->Yi, auglag->mu, auglag->Ci));
939566063dSJacob Faibussowitsch         PetscCall(VecSet(auglag->Ciwork, auglag->yi_max));
949566063dSJacob Faibussowitsch         PetscCall(VecPointwiseMin(auglag->Yi, auglag->Ciwork, auglag->Yi));
959566063dSJacob Faibussowitsch         PetscCall(VecSet(auglag->Ciwork, auglag->yi_min));
969566063dSJacob Faibussowitsch         PetscCall(VecPointwiseMax(auglag->Yi, auglag->Ciwork, auglag->Yi));
97661095bbSAlp Dener       }
98661095bbSAlp Dener       /* tolerances are updated only for non-PHR methods */
99661095bbSAlp Dener       if (auglag->type != TAO_ALMM_PHR) {
100661095bbSAlp Dener         auglag->ytol = PetscMax(tao->catol, auglag->ytol / PetscPowReal(auglag->mu, auglag->mu_pow_good));
101661095bbSAlp Dener         auglag->gtol = PetscMax(tao->gatol, auglag->gtol / auglag->mu);
102661095bbSAlp Dener       }
103661095bbSAlp Dener       updated = 1.0;
104661095bbSAlp Dener     } else {
105661095bbSAlp Dener       /* constraints are bad, update penalty factor */
106661095bbSAlp Dener       auglag->mu = PetscMin(auglag->mu_max, auglag->mu_fac * auglag->mu);
107661095bbSAlp Dener       /* tolerances are reset only for non-PHR methods */
108661095bbSAlp Dener       if (auglag->type != TAO_ALMM_PHR) {
109*c8b0c2b5SHansol Suh         auglag->ytol = PetscMax(tao->catol, 1.0 / PetscPowReal(auglag->mu, auglag->mu_pow_bad));
110661095bbSAlp Dener         auglag->gtol = PetscMax(tao->gatol, 1.0 / auglag->mu);
111661095bbSAlp Dener       }
11263a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(tao, "Penalty increased: mu = %.2g\n", (double)auglag->mu));
113661095bbSAlp Dener     }
1149566063dSJacob Faibussowitsch     PetscCall(TaoLogConvergenceHistory(tao, auglag->fval, auglag->gnorm, auglag->cnorm, tao->ksp_its));
1159566063dSJacob Faibussowitsch     PetscCall(TaoMonitor(tao, tao->niter, auglag->fval, auglag->gnorm, auglag->cnorm, updated));
116dbbe0bcdSBarry Smith     PetscUseTypeMethod(tao, convergencetest, tao->cnvP);
117661095bbSAlp Dener   }
1183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
119661095bbSAlp Dener }
120661095bbSAlp Dener 
121d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoView_ALMM(Tao tao, PetscViewer viewer)
122d71ae5a4SJacob Faibussowitsch {
123661095bbSAlp Dener   TAO_ALMM *auglag = (TAO_ALMM *)tao->data;
124661095bbSAlp Dener   PetscBool isascii;
125661095bbSAlp Dener 
126661095bbSAlp Dener   PetscFunctionBegin;
1279566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPushTab(viewer));
1289566063dSJacob Faibussowitsch   PetscCall(TaoView(auglag->subsolver, viewer));
1299566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPopTab(viewer));
1309566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
131661095bbSAlp Dener   if (isascii) {
1329566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
1339566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "ALMM Formulation Type: %s\n", TaoALMMTypes[auglag->type]));
1349566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
135661095bbSAlp Dener   }
1363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
137661095bbSAlp Dener }
138661095bbSAlp Dener 
139d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSetUp_ALMM(Tao tao)
140d71ae5a4SJacob Faibussowitsch {
141661095bbSAlp Dener   TAO_ALMM *auglag = (TAO_ALMM *)tao->data;
142661095bbSAlp Dener   VecType   vec_type;
143661095bbSAlp Dener   Vec       SL, SU;
144661095bbSAlp Dener   PetscBool is_cg = PETSC_FALSE, is_lmvm = PETSC_FALSE;
145661095bbSAlp Dener 
146661095bbSAlp Dener   PetscFunctionBegin;
14769d47153SPierre Jolivet   PetscCheck(!tao->ineq_doublesided, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "TAOALMM does not support double-sided inequality constraint definition. Please restructure your inequality constraint to fit the form c(x) >= 0.");
1483c859ba3SBarry 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.");
1499566063dSJacob Faibussowitsch   PetscCall(TaoComputeVariableBounds(tao));
150661095bbSAlp Dener   /* alias base vectors and create extras */
1519566063dSJacob Faibussowitsch   PetscCall(VecGetType(tao->solution, &vec_type));
152661095bbSAlp Dener   auglag->Px = tao->solution;
153661095bbSAlp Dener   if (!tao->gradient) { /* base gradient */
1549566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(tao->solution, &tao->gradient));
155661095bbSAlp Dener   }
156661095bbSAlp Dener   auglag->LgradX = tao->gradient;
157661095bbSAlp Dener   if (!auglag->Xwork) { /* opt var work vector */
1589566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(tao->solution, &auglag->Xwork));
159661095bbSAlp Dener   }
160661095bbSAlp Dener   if (tao->eq_constrained) {
161661095bbSAlp Dener     auglag->Ce = tao->constraints_equality;
162661095bbSAlp Dener     auglag->Ae = tao->jacobian_equality;
163661095bbSAlp Dener     if (!auglag->Ye) { /* equality multipliers */
1649566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(auglag->Ce, &auglag->Ye));
165661095bbSAlp Dener     }
16648a46eb9SPierre Jolivet     if (!auglag->Cework) PetscCall(VecDuplicate(auglag->Ce, &auglag->Cework));
167661095bbSAlp Dener   }
168661095bbSAlp Dener   if (tao->ineq_constrained) {
169661095bbSAlp Dener     auglag->Ci = tao->constraints_inequality;
170661095bbSAlp Dener     auglag->Ai = tao->jacobian_inequality;
171661095bbSAlp Dener     if (!auglag->Yi) { /* inequality multipliers */
1729566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(auglag->Ci, &auglag->Yi));
173661095bbSAlp Dener     }
17448a46eb9SPierre Jolivet     if (!auglag->Ciwork) PetscCall(VecDuplicate(auglag->Ci, &auglag->Ciwork));
175661095bbSAlp Dener     if (!auglag->Cizero) {
1769566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(auglag->Ci, &auglag->Cizero));
1779566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(auglag->Cizero));
178661095bbSAlp Dener     }
179661095bbSAlp Dener     if (!auglag->Ps) { /* slack vars */
1809566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(auglag->Ci, &auglag->Ps));
181661095bbSAlp Dener     }
182661095bbSAlp Dener     if (!auglag->LgradS) { /* slack component of Lagrangian gradient */
1839566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(auglag->Ci, &auglag->LgradS));
184661095bbSAlp Dener     }
185661095bbSAlp Dener     /* create vector for combined primal space and the associated communication objects */
186661095bbSAlp Dener     if (!auglag->P) {
1879566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(2, &auglag->Parr));
1889371c9d4SSatish Balay       auglag->Parr[0] = auglag->Px;
1899371c9d4SSatish Balay       auglag->Parr[1] = auglag->Ps;
1909566063dSJacob Faibussowitsch       PetscCall(VecConcatenate(2, auglag->Parr, &auglag->P, &auglag->Pis));
1919566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(2, &auglag->Pscatter));
1929566063dSJacob Faibussowitsch       PetscCall(VecScatterCreate(auglag->P, auglag->Pis[0], auglag->Px, NULL, &auglag->Pscatter[0]));
1939566063dSJacob Faibussowitsch       PetscCall(VecScatterCreate(auglag->P, auglag->Pis[1], auglag->Ps, NULL, &auglag->Pscatter[1]));
194661095bbSAlp Dener     }
195661095bbSAlp Dener     if (tao->eq_constrained) {
196661095bbSAlp Dener       /* create vector for combined dual space and the associated communication objects */
197661095bbSAlp Dener       if (!auglag->Y) {
1989566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(2, &auglag->Yarr));
1999371c9d4SSatish Balay         auglag->Yarr[0] = auglag->Ye;
2009371c9d4SSatish Balay         auglag->Yarr[1] = auglag->Yi;
2019566063dSJacob Faibussowitsch         PetscCall(VecConcatenate(2, auglag->Yarr, &auglag->Y, &auglag->Yis));
2029566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(2, &auglag->Yscatter));
2039566063dSJacob Faibussowitsch         PetscCall(VecScatterCreate(auglag->Y, auglag->Yis[0], auglag->Ye, NULL, &auglag->Yscatter[0]));
2049566063dSJacob Faibussowitsch         PetscCall(VecScatterCreate(auglag->Y, auglag->Yis[1], auglag->Yi, NULL, &auglag->Yscatter[1]));
205661095bbSAlp Dener       }
20648a46eb9SPierre Jolivet       if (!auglag->C) PetscCall(VecDuplicate(auglag->Y, &auglag->C));
2079371c9d4SSatish Balay     } else {
208ad540459SPierre Jolivet       if (!auglag->C) auglag->C = auglag->Ci;
209ad540459SPierre Jolivet       if (!auglag->Y) auglag->Y = auglag->Yi;
210661095bbSAlp Dener     }
211661095bbSAlp Dener   } else {
212ad540459SPierre Jolivet     if (!auglag->P) auglag->P = auglag->Px;
213ad540459SPierre Jolivet     if (!auglag->G) auglag->G = auglag->LgradX;
214ad540459SPierre Jolivet     if (!auglag->C) auglag->C = auglag->Ce;
215ad540459SPierre Jolivet     if (!auglag->Y) auglag->Y = auglag->Ye;
216661095bbSAlp Dener   }
2179a09327dSAlp Dener   /* create tao primal solution and gradient to interface with subsolver */
2189a09327dSAlp Dener   PetscCall(VecDuplicate(auglag->P, &auglag->Psub));
2199a09327dSAlp Dener   PetscCall(VecDuplicate(auglag->P, &auglag->G));
220661095bbSAlp Dener   /* initialize parameters */
221661095bbSAlp Dener   if (auglag->type == TAO_ALMM_PHR) {
222661095bbSAlp Dener     auglag->mu_fac = 10.0;
223661095bbSAlp Dener     auglag->yi_min = 0.0;
224661095bbSAlp Dener     auglag->ytol0  = 0.5;
225661095bbSAlp Dener     auglag->gtol0  = tao->gatol;
226606f75f6SBarry Smith     if (tao->gatol != tao->default_gatol && tao->catol != tao->default_catol) {
2279566063dSJacob Faibussowitsch       PetscCall(PetscInfo(tao, "TAOALMM with PHR: different gradient and constraint tolerances are not supported, setting catol = gatol\n"));
228661095bbSAlp Dener       tao->catol = tao->gatol;
229661095bbSAlp Dener     }
230661095bbSAlp Dener   }
231661095bbSAlp Dener   /* set the Lagrangian formulation type for the subsolver */
232661095bbSAlp Dener   switch (auglag->type) {
233d71ae5a4SJacob Faibussowitsch   case TAO_ALMM_CLASSIC:
234d71ae5a4SJacob Faibussowitsch     auglag->sub_obj = TaoALMMComputeAugLagAndGradient_Private;
235d71ae5a4SJacob Faibussowitsch     break;
236d71ae5a4SJacob Faibussowitsch   case TAO_ALMM_PHR:
237d71ae5a4SJacob Faibussowitsch     auglag->sub_obj = TaoALMMComputePHRLagAndGradient_Private;
238d71ae5a4SJacob Faibussowitsch     break;
239d71ae5a4SJacob Faibussowitsch   default:
240d71ae5a4SJacob Faibussowitsch     break;
241661095bbSAlp Dener   }
242661095bbSAlp Dener   /* set up the subsolver */
2439a09327dSAlp Dener   PetscCall(TaoSetSolution(auglag->subsolver, auglag->Psub));
2449566063dSJacob Faibussowitsch   PetscCall(TaoSetObjective(auglag->subsolver, TaoALMMSubsolverObjective_Private, (void *)auglag));
2459566063dSJacob Faibussowitsch   PetscCall(TaoSetObjectiveAndGradient(auglag->subsolver, NULL, TaoALMMSubsolverObjectiveAndGradient_Private, (void *)auglag));
246661095bbSAlp Dener   if (tao->bounded) {
247661095bbSAlp Dener     /* make sure that the subsolver is a bound-constrained method */
2489566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)auglag->subsolver, TAOCG, &is_cg));
2499566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)auglag->subsolver, TAOLMVM, &is_lmvm));
250661095bbSAlp Dener     if (is_cg) {
2519566063dSJacob Faibussowitsch       PetscCall(TaoSetType(auglag->subsolver, TAOBNCG));
2529d3446b2SPierre Jolivet       PetscCall(PetscInfo(tao, "TAOCG detected for bound-constrained problem, switching to TAOBNCG.\n"));
253661095bbSAlp Dener     }
254661095bbSAlp Dener     if (is_lmvm) {
2559566063dSJacob Faibussowitsch       PetscCall(TaoSetType(auglag->subsolver, TAOBQNLS));
2569d3446b2SPierre Jolivet       PetscCall(PetscInfo(tao, "TAOLMVM detected for bound-constrained problem, switching to TAOBQNLS.\n"));
257661095bbSAlp Dener     }
258661095bbSAlp Dener     /* create lower and upper bound clone vectors for subsolver */
25948a46eb9SPierre Jolivet     if (!auglag->PL) PetscCall(VecDuplicate(auglag->P, &auglag->PL));
26048a46eb9SPierre Jolivet     if (!auglag->PU) PetscCall(VecDuplicate(auglag->P, &auglag->PU));
261661095bbSAlp Dener     if (tao->ineq_constrained) {
262661095bbSAlp Dener       /* create lower and upper bounds for slack, set lower to 0 */
2639566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(auglag->Ci, &SL));
2649566063dSJacob Faibussowitsch       PetscCall(VecSet(SL, 0.0));
2659566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(auglag->Ci, &SU));
2669566063dSJacob Faibussowitsch       PetscCall(VecSet(SU, PETSC_INFINITY));
267661095bbSAlp Dener       /* combine opt var bounds with slack bounds */
2689566063dSJacob Faibussowitsch       PetscCall(TaoALMMCombinePrimal_Private(tao, tao->XL, SL, auglag->PL));
2699566063dSJacob Faibussowitsch       PetscCall(TaoALMMCombinePrimal_Private(tao, tao->XU, SU, auglag->PU));
270661095bbSAlp Dener       /* destroy work vectors */
2719566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&SL));
2729566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&SU));
273661095bbSAlp Dener     } else {
274661095bbSAlp Dener       /* no inequality constraints, just copy bounds into the subsolver */
2759566063dSJacob Faibussowitsch       PetscCall(VecCopy(tao->XL, auglag->PL));
2769566063dSJacob Faibussowitsch       PetscCall(VecCopy(tao->XU, auglag->PU));
277661095bbSAlp Dener     }
2789566063dSJacob Faibussowitsch     PetscCall(TaoSetVariableBounds(auglag->subsolver, auglag->PL, auglag->PU));
279*c8b0c2b5SHansol Suh   } else {
280*c8b0c2b5SHansol Suh     /* CLASSIC's slack variable is bounded, so need to set bounds */
281*c8b0c2b5SHansol Suh     //TODO what happens for non-constrained ALMM CLASSIC?
282*c8b0c2b5SHansol Suh     if (auglag->type == TAO_ALMM_CLASSIC) {
283*c8b0c2b5SHansol Suh       if (tao->ineq_constrained) {
284*c8b0c2b5SHansol Suh         /* create lower and upper bound clone vectors for subsolver
285*c8b0c2b5SHansol Suh          * They should be NFINITY and INFINITY                       */
286*c8b0c2b5SHansol Suh         if (!auglag->PL) PetscCall(VecDuplicate(auglag->P, &auglag->PL));
287*c8b0c2b5SHansol Suh         if (!auglag->PU) PetscCall(VecDuplicate(auglag->P, &auglag->PU));
288*c8b0c2b5SHansol Suh         PetscCall(VecSet(auglag->PL, PETSC_NINFINITY));
289*c8b0c2b5SHansol Suh         PetscCall(VecSet(auglag->PU, PETSC_INFINITY));
290*c8b0c2b5SHansol Suh         /* create lower and upper bounds for slack, set lower to 0 */
291*c8b0c2b5SHansol Suh         PetscCall(VecDuplicate(auglag->Ci, &SL));
292*c8b0c2b5SHansol Suh         PetscCall(VecSet(SL, 0.0));
293*c8b0c2b5SHansol Suh         PetscCall(VecDuplicate(auglag->Ci, &SU));
294*c8b0c2b5SHansol Suh         PetscCall(VecSet(SU, PETSC_INFINITY));
295*c8b0c2b5SHansol Suh         /* PL, PU is already set. Only copy Slack variable parts */
296*c8b0c2b5SHansol Suh         PetscCall(VecScatterBegin(auglag->Pscatter[1], SL, auglag->PL, INSERT_VALUES, SCATTER_REVERSE));
297*c8b0c2b5SHansol Suh         PetscCall(VecScatterEnd(auglag->Pscatter[1], SL, auglag->PL, INSERT_VALUES, SCATTER_REVERSE));
298*c8b0c2b5SHansol Suh         PetscCall(VecScatterBegin(auglag->Pscatter[1], SU, auglag->PU, INSERT_VALUES, SCATTER_REVERSE));
299*c8b0c2b5SHansol Suh         PetscCall(VecScatterEnd(auglag->Pscatter[1], SU, auglag->PU, INSERT_VALUES, SCATTER_REVERSE));
300*c8b0c2b5SHansol Suh         /* destroy work vectors */
301*c8b0c2b5SHansol Suh         PetscCall(VecDestroy(&SL));
302*c8b0c2b5SHansol Suh         PetscCall(VecDestroy(&SU));
303*c8b0c2b5SHansol Suh         /* make sure that the subsolver is a bound-constrained method
304*c8b0c2b5SHansol Suh          * Unfortunately duplicate code                                 */
305*c8b0c2b5SHansol Suh         PetscCall(PetscObjectTypeCompare((PetscObject)auglag->subsolver, TAOCG, &is_cg));
306*c8b0c2b5SHansol Suh         PetscCall(PetscObjectTypeCompare((PetscObject)auglag->subsolver, TAOLMVM, &is_lmvm));
307*c8b0c2b5SHansol Suh         if (is_cg) {
308*c8b0c2b5SHansol Suh           PetscCall(TaoSetType(auglag->subsolver, TAOBNCG));
309*c8b0c2b5SHansol Suh           PetscCall(PetscInfo(tao, "TAOCG detected for TAO_ALMM_CLASSIC, switching to TAOBNCG.\n"));
310*c8b0c2b5SHansol Suh         }
311*c8b0c2b5SHansol Suh         if (is_lmvm) {
312*c8b0c2b5SHansol Suh           PetscCall(TaoSetType(auglag->subsolver, TAOBQNLS));
313*c8b0c2b5SHansol Suh           PetscCall(PetscInfo(tao, "TAOLMVM detected for TAO_ALMM_CLASSIC, switching to TAOBQNLS.\n"));
314*c8b0c2b5SHansol Suh         }
315*c8b0c2b5SHansol Suh         PetscCall(TaoSetVariableBounds(auglag->subsolver, auglag->PL, auglag->PU));
316*c8b0c2b5SHansol Suh       }
317*c8b0c2b5SHansol Suh     }
318661095bbSAlp Dener   }
3199566063dSJacob Faibussowitsch   PetscCall(TaoSetUp(auglag->subsolver));
3203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
321661095bbSAlp Dener }
322661095bbSAlp Dener 
323d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoDestroy_ALMM(Tao tao)
324d71ae5a4SJacob Faibussowitsch {
325661095bbSAlp Dener   TAO_ALMM *auglag = (TAO_ALMM *)tao->data;
326661095bbSAlp Dener 
327661095bbSAlp Dener   PetscFunctionBegin;
3289566063dSJacob Faibussowitsch   PetscCall(TaoDestroy(&auglag->subsolver));
3299a09327dSAlp Dener   PetscCall(VecDestroy(&auglag->Psub));
3309a09327dSAlp Dener   PetscCall(VecDestroy(&auglag->G));
331661095bbSAlp Dener   if (tao->setupcalled) {
3329a09327dSAlp Dener     PetscCall(VecDestroy(&auglag->Xwork));
333661095bbSAlp Dener     if (tao->eq_constrained) {
3349566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&auglag->Ye));     /* equality multipliers */
3359566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&auglag->Cework)); /* equality work vector */
336661095bbSAlp Dener     }
337661095bbSAlp Dener     if (tao->ineq_constrained) {
3389566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&auglag->Ps));     /* slack vars */
3399566063dSJacob Faibussowitsch       PetscCall(PetscFree(auglag->Parr));     /* array of primal vectors */
3409566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&auglag->LgradS)); /* slack grad */
3419566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&auglag->Cizero)); /* zero vector for pointwise max */
3429566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&auglag->Yi));     /* inequality multipliers */
3439566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&auglag->Ciwork)); /* inequality work vector */
3449a09327dSAlp Dener       PetscCall(VecDestroy(&auglag->P));      /* combo primal solution */
3459566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&auglag->Pis[0]));  /* index set for X inside P */
3469566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&auglag->Pis[1]));  /* index set for S inside P */
3479566063dSJacob Faibussowitsch       PetscCall(PetscFree(auglag->Pis));      /* array of P index sets */
3489566063dSJacob Faibussowitsch       PetscCall(VecScatterDestroy(&auglag->Pscatter[0]));
3499566063dSJacob Faibussowitsch       PetscCall(VecScatterDestroy(&auglag->Pscatter[1]));
3509566063dSJacob Faibussowitsch       PetscCall(PetscFree(auglag->Pscatter));
351661095bbSAlp Dener       if (tao->eq_constrained) {
3529566063dSJacob Faibussowitsch         PetscCall(VecDestroy(&auglag->Y));     /* combo multipliers */
3539566063dSJacob Faibussowitsch         PetscCall(PetscFree(auglag->Yarr));    /* array of dual vectors */
3549566063dSJacob Faibussowitsch         PetscCall(VecDestroy(&auglag->C));     /* combo constraints */
3559566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&auglag->Yis[0])); /* index set for Ye inside Y */
3569566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&auglag->Yis[1])); /* index set for Yi inside Y */
3579566063dSJacob Faibussowitsch         PetscCall(PetscFree(auglag->Yis));
3589566063dSJacob Faibussowitsch         PetscCall(VecScatterDestroy(&auglag->Yscatter[0]));
3599566063dSJacob Faibussowitsch         PetscCall(VecScatterDestroy(&auglag->Yscatter[1]));
3609566063dSJacob Faibussowitsch         PetscCall(PetscFree(auglag->Yscatter));
361661095bbSAlp Dener       }
362661095bbSAlp Dener     }
363661095bbSAlp Dener     if (tao->bounded) {
3649566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&auglag->PL)); /* lower bounds for subsolver */
3659566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&auglag->PU)); /* upper bounds for subsolver */
366*c8b0c2b5SHansol Suh     } else {
367*c8b0c2b5SHansol Suh       if (auglag->type == TAO_ALMM_CLASSIC) {
368*c8b0c2b5SHansol Suh         PetscCall(VecDestroy(&auglag->PL)); /* lower bounds for subsolver */
369*c8b0c2b5SHansol Suh         PetscCall(VecDestroy(&auglag->PU)); /* upper bounds for subsolver */
370*c8b0c2b5SHansol Suh       }
371661095bbSAlp Dener     }
372661095bbSAlp Dener   }
3732e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetType_C", NULL));
3742e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMSetType_C", NULL));
3752e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetSubsolver_C", NULL));
3762e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMSetSubsolver_C", NULL));
3772e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetMultipliers_C", NULL));
3782e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMSetMultipliers_C", NULL));
3792e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetPrimalIS_C", NULL));
3802e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetDualIS_C", NULL));
3819566063dSJacob Faibussowitsch   PetscCall(PetscFree(tao->data));
3823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
383661095bbSAlp Dener }
384661095bbSAlp Dener 
385ce78bad3SBarry Smith static PetscErrorCode TaoSetFromOptions_ALMM(Tao tao, PetscOptionItems PetscOptionsObject)
386d71ae5a4SJacob Faibussowitsch {
387661095bbSAlp Dener   TAO_ALMM *auglag = (TAO_ALMM *)tao->data;
388661095bbSAlp Dener   PetscInt  i;
389661095bbSAlp Dener 
390661095bbSAlp Dener   PetscFunctionBegin;
39171075aafSSatish Balay   PetscOptionsHeadBegin(PetscOptionsObject, "Augmented Lagrangian multiplier method solves problems with general constraints by converting them into a sequence of unconstrained problems.");
3929566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_almm_mu_init", "initial penalty parameter", "", auglag->mu0, &auglag->mu0, NULL));
3939566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_almm_mu_factor", "increase factor for the penalty parameter", "", auglag->mu_fac, &auglag->mu_fac, NULL));
3949566063dSJacob 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));
3959566063dSJacob 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));
3969566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_almm_mu_max", "maximum safeguard for penalty parameter updates", "", auglag->mu_max, &auglag->mu_max, NULL));
3979566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_almm_ye_min", "minimum safeguard for equality multiplier updates", "", auglag->ye_min, &auglag->ye_min, NULL));
3989566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_almm_ye_max", "maximum safeguard for equality multipliers updates", "", auglag->ye_max, &auglag->ye_max, NULL));
3999566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_almm_yi_min", "minimum safeguard for inequality multipliers updates", "", auglag->yi_min, &auglag->yi_min, NULL));
4009566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_almm_yi_max", "maximum safeguard for inequality multipliers updates", "", auglag->yi_max, &auglag->yi_max, NULL));
4019566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-tao_almm_type", "augmented Lagrangian formulation type for the subproblem", "TaoALMMType", TaoALMMTypes, (PetscEnum)auglag->type, (PetscEnum *)&auglag->type, NULL));
402d0609cedSBarry Smith   PetscOptionsHeadEnd();
4039566063dSJacob Faibussowitsch   PetscCall(TaoSetOptionsPrefix(auglag->subsolver, ((PetscObject)tao)->prefix));
4049566063dSJacob Faibussowitsch   PetscCall(TaoAppendOptionsPrefix(auglag->subsolver, "tao_almm_subsolver_"));
4059566063dSJacob Faibussowitsch   PetscCall(TaoSetFromOptions(auglag->subsolver));
406661095bbSAlp Dener   for (i = 0; i < tao->numbermonitors; i++) {
4079566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)tao->monitorcontext[i]));
40810978b7dSBarry Smith     PetscCall(TaoMonitorSet(auglag->subsolver, tao->monitor[i], tao->monitorcontext[i], tao->monitordestroy[i]));
40910978b7dSBarry Smith     if (tao->monitor[i] == TaoMonitorDefault || tao->monitor[i] == TaoMonitorConstraintNorm || tao->monitor[i] == TaoMonitorGlobalization || tao->monitor[i] == TaoMonitorDefaultShort) auglag->info = PETSC_TRUE;
410661095bbSAlp Dener   }
4113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
412661095bbSAlp Dener }
413661095bbSAlp Dener 
414661095bbSAlp Dener /* -------------------------------------------------------- */
415661095bbSAlp Dener 
416661095bbSAlp Dener /*MC
417c0298470SBarry Smith   TAOALMM - Augmented Lagrangian multiplier method for solving nonlinear optimization problems with general constraints.
418661095bbSAlp Dener 
419661095bbSAlp Dener   Options Database Keys:
420661095bbSAlp Dener + -tao_almm_mu_init <real>       - initial penalty parameter (default: 10.)
421661095bbSAlp Dener . -tao_almm_mu_factor <real>     - increase factor for the penalty parameter (default: 100.)
422661095bbSAlp Dener . -tao_almm_mu_max <real>        - maximum safeguard for penalty parameter updates (default: 1.e20)
423661095bbSAlp Dener . -tao_almm_mu_power_good <real> - exponential for penalty parameter when multiplier update is accepted (default: 0.9)
424661095bbSAlp Dener . -tao_almm_mu_power_bad <real>  - exponential for penalty parameter when multiplier update is rejected (default: 0.1)
425661095bbSAlp Dener . -tao_almm_ye_min <real>        - minimum safeguard for equality multiplier updates (default: -1.e20)
426661095bbSAlp Dener . -tao_almm_ye_max <real>        - maximum safeguard for equality multiplier updates (default: 1.e20)
427661095bbSAlp Dener . -tao_almm_yi_min <real>        - minimum safeguard for inequality multiplier updates (default: -1.e20)
428661095bbSAlp Dener . -tao_almm_yi_max <real>        - maximum safeguard for inequality multiplier updates (default: 1.e20)
4299a09327dSAlp Dener - -tao_almm_type <phr,classic>   - change formulation of the augmented Lagrangian merit function for the subproblem (default: phr)
430661095bbSAlp Dener 
431661095bbSAlp Dener   Level: beginner
432661095bbSAlp Dener 
433661095bbSAlp Dener   Notes:
434661095bbSAlp Dener   This method converts a constrained problem into a sequence of unconstrained problems via the augmented
435661095bbSAlp Dener   Lagrangian merit function. Bound constraints are pushed down to the subproblem without any modifications.
436661095bbSAlp Dener 
437661095bbSAlp Dener   Two formulations are offered for the subproblem: canonical Hestenes-Powell augmented Lagrangian with slack
438661095bbSAlp Dener   variables for inequality constraints, and a slack-less Powell-Hestenes-Rockafellar (PHR) formulation utilizing a
4399a09327dSAlp Dener   pointwise max() penalty on inequality constraints. The canonical augmented Lagrangian formulation may converge
4409a09327dSAlp Dener   faster for smaller problems but is highly susceptible to poor step lengths in the subproblem due to the positivity
4419a09327dSAlp Dener   constraint on slack variables. PHR avoids this issue by eliminating the slack variables entirely, and is highly
4429a09327dSAlp Dener   desirable for problems with a large number of inequality constraints.
443661095bbSAlp Dener 
4449a09327dSAlp Dener   The subproblem is solved using a nested first-order TAO solver (default: `TAOBQNLS`). The user can retrieve a
4459a09327dSAlp Dener   pointer to the subsolver via `TaoALMMGetSubsolver()` or pass command line arguments to it using the
4469a09327dSAlp Dener   "-tao_almm_subsolver_" prefix. Currently, `TAOALMM` does not support second-order methods for the subproblem.
447661095bbSAlp Dener 
448661095bbSAlp Dener .vb
449661095bbSAlp Dener   while unconverged
450661095bbSAlp Dener     solve argmin_x L(x) s.t. l <= x <= u
451661095bbSAlp Dener     if ||c|| <= y_tol
452661095bbSAlp Dener       if ||c|| <= c_tol && ||Lgrad|| <= g_tol:
453661095bbSAlp Dener         problem converged, return solution
454661095bbSAlp Dener       else
455661095bbSAlp Dener         constraints sufficiently improved
456661095bbSAlp Dener         update multipliers and tighten tolerances
457661095bbSAlp Dener       endif
458661095bbSAlp Dener     else
459661095bbSAlp Dener       constraints did not improve
460661095bbSAlp Dener       update penalty and loosen tolerances
461661095bbSAlp Dener     endif
462661095bbSAlp Dener   endwhile
463661095bbSAlp Dener .ve
464661095bbSAlp Dener 
465c0298470SBarry Smith .seealso: `TAOALMM`, `Tao`, `TaoALMMGetType()`, `TaoALMMSetType()`, `TaoALMMSetSubsolver()`, `TaoALMMGetSubsolver()`,
466db781477SPatrick Sanan           `TaoALMMGetMultipliers()`, `TaoALMMSetMultipliers()`, `TaoALMMGetPrimalIS()`, `TaoALMMGetDualIS()`
467661095bbSAlp Dener M*/
468d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TaoCreate_ALMM(Tao tao)
469d71ae5a4SJacob Faibussowitsch {
470661095bbSAlp Dener   TAO_ALMM *auglag;
471661095bbSAlp Dener 
472661095bbSAlp Dener   PetscFunctionBegin;
4734dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&auglag));
474661095bbSAlp Dener 
475661095bbSAlp Dener   tao->ops->destroy        = TaoDestroy_ALMM;
476661095bbSAlp Dener   tao->ops->setup          = TaoSetUp_ALMM;
477661095bbSAlp Dener   tao->ops->setfromoptions = TaoSetFromOptions_ALMM;
478661095bbSAlp Dener   tao->ops->view           = TaoView_ALMM;
479661095bbSAlp Dener   tao->ops->solve          = TaoSolve_ALMM;
480661095bbSAlp Dener 
481606f75f6SBarry Smith   PetscCall(TaoParametersInitialize(tao));
482606f75f6SBarry Smith   PetscObjectParameterSetDefault(tao, gatol, 1.e-5);
483606f75f6SBarry Smith   PetscObjectParameterSetDefault(tao, grtol, 0.0);
484606f75f6SBarry Smith   PetscObjectParameterSetDefault(tao, gttol, 0.0);
485606f75f6SBarry Smith   PetscObjectParameterSetDefault(tao, catol, 1.e-5);
486606f75f6SBarry Smith   PetscObjectParameterSetDefault(tao, crtol, 0.0);
487661095bbSAlp Dener 
488661095bbSAlp Dener   tao->data           = (void *)auglag;
489661095bbSAlp Dener   auglag->parent      = tao;
490661095bbSAlp Dener   auglag->mu0         = 10.0;
491661095bbSAlp Dener   auglag->mu          = auglag->mu0;
492661095bbSAlp Dener   auglag->mu_fac      = 10.0;
493661095bbSAlp Dener   auglag->mu_max      = PETSC_INFINITY;
494661095bbSAlp Dener   auglag->mu_pow_good = 0.9;
495661095bbSAlp Dener   auglag->mu_pow_bad  = 0.1;
496661095bbSAlp Dener   auglag->ye_min      = PETSC_NINFINITY;
497661095bbSAlp Dener   auglag->ye_max      = PETSC_INFINITY;
498661095bbSAlp Dener   auglag->yi_min      = PETSC_NINFINITY;
499661095bbSAlp Dener   auglag->yi_max      = PETSC_INFINITY;
500*c8b0c2b5SHansol Suh   auglag->ytol0       = 1.0 / PetscPowReal(auglag->mu0, auglag->mu_pow_bad);
501661095bbSAlp Dener   auglag->ytol        = auglag->ytol0;
502661095bbSAlp Dener   auglag->gtol0       = 1.0 / auglag->mu0;
503661095bbSAlp Dener   auglag->gtol        = auglag->gtol0;
504661095bbSAlp Dener 
505661095bbSAlp Dener   auglag->sub_obj = TaoALMMComputeAugLagAndGradient_Private;
5069a09327dSAlp Dener   auglag->type    = TAO_ALMM_PHR;
507661095bbSAlp Dener   auglag->info    = PETSC_FALSE;
508661095bbSAlp Dener 
5099566063dSJacob Faibussowitsch   PetscCall(TaoCreate(PetscObjectComm((PetscObject)tao), &auglag->subsolver));
5109a09327dSAlp Dener   PetscCall(TaoSetType(auglag->subsolver, TAOBQNLS));
5119566063dSJacob Faibussowitsch   PetscCall(TaoSetTolerances(auglag->subsolver, auglag->gtol, 0.0, 0.0));
5129566063dSJacob Faibussowitsch   PetscCall(TaoSetMaximumIterations(auglag->subsolver, 1000));
5139566063dSJacob Faibussowitsch   PetscCall(TaoSetMaximumFunctionEvaluations(auglag->subsolver, 10000));
5149566063dSJacob Faibussowitsch   PetscCall(TaoSetFunctionLowerBound(auglag->subsolver, PETSC_NINFINITY));
5159566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)auglag->subsolver, (PetscObject)tao, 1));
516661095bbSAlp Dener 
5179566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetType_C", TaoALMMGetType_Private));
5189566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMSetType_C", TaoALMMSetType_Private));
5199566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetSubsolver_C", TaoALMMGetSubsolver_Private));
5209566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMSetSubsolver_C", TaoALMMSetSubsolver_Private));
5219566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetMultipliers_C", TaoALMMGetMultipliers_Private));
5229566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMSetMultipliers_C", TaoALMMSetMultipliers_Private));
5239566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetPrimalIS_C", TaoALMMGetPrimalIS_Private));
5249566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetDualIS_C", TaoALMMGetDualIS_Private));
5253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
526661095bbSAlp Dener }
527661095bbSAlp Dener 
528d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoALMMCombinePrimal_Private(Tao tao, Vec X, Vec S, Vec P)
529d71ae5a4SJacob Faibussowitsch {
530661095bbSAlp Dener   TAO_ALMM *auglag = (TAO_ALMM *)tao->data;
531661095bbSAlp Dener 
532661095bbSAlp Dener   PetscFunctionBegin;
533661095bbSAlp Dener   if (tao->ineq_constrained) {
5349566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(auglag->Pscatter[0], X, P, INSERT_VALUES, SCATTER_REVERSE));
5359566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(auglag->Pscatter[0], X, P, INSERT_VALUES, SCATTER_REVERSE));
5369566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(auglag->Pscatter[1], S, P, INSERT_VALUES, SCATTER_REVERSE));
5379566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(auglag->Pscatter[1], S, P, INSERT_VALUES, SCATTER_REVERSE));
538661095bbSAlp Dener   } else {
5399566063dSJacob Faibussowitsch     PetscCall(VecCopy(X, P));
540661095bbSAlp Dener   }
5413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
542661095bbSAlp Dener }
543661095bbSAlp Dener 
544d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoALMMCombineDual_Private(Tao tao, Vec EQ, Vec IN, Vec Y)
545d71ae5a4SJacob Faibussowitsch {
546661095bbSAlp Dener   TAO_ALMM *auglag = (TAO_ALMM *)tao->data;
547661095bbSAlp Dener 
548661095bbSAlp Dener   PetscFunctionBegin;
549661095bbSAlp Dener   if (tao->eq_constrained) {
550661095bbSAlp Dener     if (tao->ineq_constrained) {
5519566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(auglag->Yscatter[0], EQ, Y, INSERT_VALUES, SCATTER_REVERSE));
5529566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(auglag->Yscatter[0], EQ, Y, INSERT_VALUES, SCATTER_REVERSE));
5539566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(auglag->Yscatter[1], IN, Y, INSERT_VALUES, SCATTER_REVERSE));
5549566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(auglag->Yscatter[1], IN, Y, INSERT_VALUES, SCATTER_REVERSE));
555661095bbSAlp Dener     } else {
5569566063dSJacob Faibussowitsch       PetscCall(VecCopy(EQ, Y));
557661095bbSAlp Dener     }
558661095bbSAlp Dener   } else {
5599566063dSJacob Faibussowitsch     PetscCall(VecCopy(IN, Y));
560661095bbSAlp Dener   }
5613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
562661095bbSAlp Dener }
563661095bbSAlp Dener 
564d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoALMMSplitPrimal_Private(Tao tao, Vec P, Vec X, Vec S)
565d71ae5a4SJacob Faibussowitsch {
566661095bbSAlp Dener   TAO_ALMM *auglag = (TAO_ALMM *)tao->data;
567661095bbSAlp Dener 
568661095bbSAlp Dener   PetscFunctionBegin;
569661095bbSAlp Dener   if (tao->ineq_constrained) {
5709566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(auglag->Pscatter[0], P, X, INSERT_VALUES, SCATTER_FORWARD));
5719566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(auglag->Pscatter[0], P, X, INSERT_VALUES, SCATTER_FORWARD));
5729566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(auglag->Pscatter[1], P, S, INSERT_VALUES, SCATTER_FORWARD));
5739566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(auglag->Pscatter[1], P, S, INSERT_VALUES, SCATTER_FORWARD));
574661095bbSAlp Dener   } else {
5759566063dSJacob Faibussowitsch     PetscCall(VecCopy(P, X));
576661095bbSAlp Dener   }
5773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
578661095bbSAlp Dener }
579661095bbSAlp Dener 
580661095bbSAlp Dener /* this assumes that the latest constraints are stored in Ce and Ci, and also combined in C */
581d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoALMMComputeOptimalityNorms_Private(Tao tao)
582d71ae5a4SJacob Faibussowitsch {
583661095bbSAlp Dener   TAO_ALMM *auglag = (TAO_ALMM *)tao->data;
584661095bbSAlp Dener 
585661095bbSAlp Dener   PetscFunctionBegin;
586661095bbSAlp Dener   /* if bounded, project the gradient */
5871baa6e33SBarry Smith   if (tao->bounded) PetscCall(VecBoundGradientProjection(auglag->LgradX, auglag->Px, tao->XL, tao->XU, auglag->LgradX));
588661095bbSAlp Dener   if (auglag->type == TAO_ALMM_PHR) {
5899566063dSJacob Faibussowitsch     PetscCall(VecNorm(auglag->LgradX, NORM_INFINITY, &auglag->gnorm));
590661095bbSAlp Dener     auglag->cenorm = 0.0;
59148a46eb9SPierre Jolivet     if (tao->eq_constrained) PetscCall(VecNorm(auglag->Ce, NORM_INFINITY, &auglag->cenorm));
592661095bbSAlp Dener     auglag->cinorm = 0.0;
593661095bbSAlp Dener     if (tao->ineq_constrained) {
5949566063dSJacob Faibussowitsch       PetscCall(VecCopy(auglag->Yi, auglag->Ciwork));
5959566063dSJacob Faibussowitsch       PetscCall(VecScale(auglag->Ciwork, -1.0 / auglag->mu));
5969566063dSJacob Faibussowitsch       PetscCall(VecPointwiseMax(auglag->Ciwork, auglag->Ci, auglag->Ciwork));
5979566063dSJacob Faibussowitsch       PetscCall(VecNorm(auglag->Ciwork, NORM_INFINITY, &auglag->cinorm));
598661095bbSAlp Dener     }
599661095bbSAlp Dener     auglag->cnorm_old = auglag->cnorm;
600661095bbSAlp Dener     auglag->cnorm     = PetscMax(auglag->cenorm, auglag->cinorm);
601661095bbSAlp Dener     auglag->ytol      = auglag->ytol0 * auglag->cnorm_old;
602661095bbSAlp Dener   } else {
6039566063dSJacob Faibussowitsch     PetscCall(VecNorm(auglag->LgradX, NORM_2, &auglag->gnorm));
6049566063dSJacob Faibussowitsch     PetscCall(VecNorm(auglag->C, NORM_2, &auglag->cnorm));
605661095bbSAlp Dener   }
6063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
607661095bbSAlp Dener }
608661095bbSAlp Dener 
6099a09327dSAlp Dener static PetscErrorCode TaoALMMEvaluateIterate_Private(Tao tao)
610d71ae5a4SJacob Faibussowitsch {
611661095bbSAlp Dener   TAO_ALMM *auglag = (TAO_ALMM *)tao->data;
612661095bbSAlp Dener 
613661095bbSAlp Dener   PetscFunctionBegin;
614661095bbSAlp Dener   /* split solution into primal and slack components */
6159566063dSJacob Faibussowitsch   PetscCall(TaoALMMSplitPrimal_Private(tao, auglag->P, auglag->Px, auglag->Ps));
616661095bbSAlp Dener 
617661095bbSAlp Dener   /* compute f, df/dx and the constraints */
6189566063dSJacob Faibussowitsch   PetscCall(TaoComputeObjectiveAndGradient(tao, auglag->Px, &auglag->fval, auglag->LgradX));
619661095bbSAlp Dener   if (tao->eq_constrained) {
6209566063dSJacob Faibussowitsch     PetscCall(TaoComputeEqualityConstraints(tao, auglag->Px, auglag->Ce));
6219566063dSJacob Faibussowitsch     PetscCall(TaoComputeJacobianEquality(tao, auglag->Px, auglag->Ae, auglag->Ae));
622661095bbSAlp Dener   }
623661095bbSAlp Dener   if (tao->ineq_constrained) {
6249566063dSJacob Faibussowitsch     PetscCall(TaoComputeInequalityConstraints(tao, auglag->Px, auglag->Ci));
6259566063dSJacob Faibussowitsch     PetscCall(TaoComputeJacobianInequality(tao, auglag->Px, auglag->Ai, auglag->Ai));
626661095bbSAlp Dener     switch (auglag->type) {
6277721a36fSStefano Zampini     case TAO_ALMM_CLASSIC:
628661095bbSAlp Dener       /* classic formulation converts inequality to equality constraints via slack variables */
6299566063dSJacob Faibussowitsch       PetscCall(VecAXPY(auglag->Ci, -1.0, auglag->Ps));
630661095bbSAlp Dener       break;
6317721a36fSStefano Zampini     case TAO_ALMM_PHR:
632661095bbSAlp Dener       /* PHR is based on Ci <= 0 while TAO defines Ci >= 0 so we hit it with a negative sign */
6339566063dSJacob Faibussowitsch       PetscCall(VecScale(auglag->Ci, -1.0));
6349566063dSJacob Faibussowitsch       PetscCall(MatScale(auglag->Ai, -1.0));
635661095bbSAlp Dener       break;
636d71ae5a4SJacob Faibussowitsch     default:
637d71ae5a4SJacob Faibussowitsch       break;
638661095bbSAlp Dener     }
639661095bbSAlp Dener   }
640661095bbSAlp Dener   /* combine constraints into one vector */
6419566063dSJacob Faibussowitsch   PetscCall(TaoALMMCombineDual_Private(tao, auglag->Ce, auglag->Ci, auglag->C));
6423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
643661095bbSAlp Dener }
644661095bbSAlp Dener 
645661095bbSAlp Dener /*
646e2f36375SPaul T. Kühner Lphr = f + 0.5*mu*[ (Ce + Ye/mu)^T (Ce + Ye/mu) + pmax(0, Ci + Yi/mu)^T pmax(0, Ci + Yi/mu)]
647661095bbSAlp Dener 
648e2f36375SPaul T. Kühner dLphr/dX = dF/dX + mu*[ (Ce + Ye/mu)^T Ae + pmax(0, Ci + Yi/mu)^T Ai]
649661095bbSAlp Dener 
650661095bbSAlp Dener dLphr/dS = 0
651661095bbSAlp Dener */
652d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoALMMComputePHRLagAndGradient_Private(Tao tao)
653d71ae5a4SJacob Faibussowitsch {
654661095bbSAlp Dener   TAO_ALMM *auglag  = (TAO_ALMM *)tao->data;
655661095bbSAlp Dener   PetscReal eq_norm = 0.0, ineq_norm = 0.0;
656661095bbSAlp Dener 
657661095bbSAlp Dener   PetscFunctionBegin;
6589a09327dSAlp Dener   PetscCall(TaoALMMEvaluateIterate_Private(tao));
659661095bbSAlp Dener   if (tao->eq_constrained) {
660661095bbSAlp Dener     /* Ce_work = mu*(Ce + Ye/mu) */
6619566063dSJacob Faibussowitsch     PetscCall(VecWAXPY(auglag->Cework, 1.0 / auglag->mu, auglag->Ye, auglag->Ce));
6629566063dSJacob Faibussowitsch     PetscCall(VecDot(auglag->Cework, auglag->Cework, &eq_norm)); /* contribution to scalar Lagrangian */
6639566063dSJacob Faibussowitsch     PetscCall(VecScale(auglag->Cework, auglag->mu));
664661095bbSAlp Dener     /* dL/dX += mu*(Ce + Ye/mu)^T Ae */
6659566063dSJacob Faibussowitsch     PetscCall(MatMultTransposeAdd(auglag->Ae, auglag->Cework, auglag->LgradX, auglag->LgradX));
666661095bbSAlp Dener   }
667661095bbSAlp Dener   if (tao->ineq_constrained) {
668661095bbSAlp Dener     /* Ci_work = mu * pmax(0, Ci + Yi/mu) where pmax() is pointwise max() */
6699566063dSJacob Faibussowitsch     PetscCall(VecWAXPY(auglag->Ciwork, 1.0 / auglag->mu, auglag->Yi, auglag->Ci));
6709566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMax(auglag->Ciwork, auglag->Cizero, auglag->Ciwork));
6719566063dSJacob Faibussowitsch     PetscCall(VecDot(auglag->Ciwork, auglag->Ciwork, &ineq_norm)); /* contribution to scalar Lagrangian */
672661095bbSAlp Dener     /* dL/dX += mu * pmax(0, Ci + Yi/mu)^T Ai */
6739566063dSJacob Faibussowitsch     PetscCall(VecScale(auglag->Ciwork, auglag->mu));
6749566063dSJacob Faibussowitsch     PetscCall(MatMultTransposeAdd(auglag->Ai, auglag->Ciwork, auglag->LgradX, auglag->LgradX));
675a5b23f4aSJose E. Roman     /* dL/dS = 0 because there are no slacks in PHR */
6769566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(auglag->LgradS));
677661095bbSAlp Dener   }
678661095bbSAlp Dener   /* combine gradient together */
6799566063dSJacob Faibussowitsch   PetscCall(TaoALMMCombinePrimal_Private(tao, auglag->LgradX, auglag->LgradS, auglag->G));
680661095bbSAlp 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)] */
6815f80ce2aSJacob Faibussowitsch   auglag->Lval = auglag->fval + 0.5 * auglag->mu * (eq_norm + ineq_norm);
6823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
683661095bbSAlp Dener }
684661095bbSAlp Dener 
685661095bbSAlp Dener /*
686661095bbSAlp Dener Lc = F + Ye^TCe + Yi^T(Ci - S) + 0.5*mu*[Ce^TCe + (Ci - S)^T(Ci - S)]
687661095bbSAlp Dener 
688a9b58be2SPaul T. Kühner dLc/dX = dF/dX + Ye^TAe + Yi^TAi + mu*[Ce^TAe + (Ci - S)^TAi]
689661095bbSAlp Dener 
690661095bbSAlp Dener dLc/dS = -[Yi + mu*(Ci - S)]
691661095bbSAlp Dener */
692d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoALMMComputeAugLagAndGradient_Private(Tao tao)
693d71ae5a4SJacob Faibussowitsch {
694661095bbSAlp Dener   TAO_ALMM *auglag = (TAO_ALMM *)tao->data;
695661095bbSAlp Dener   PetscReal yeTce = 0.0, yiTcims = 0.0, ceTce = 0.0, cimsTcims = 0.0;
696661095bbSAlp Dener 
697661095bbSAlp Dener   PetscFunctionBegin;
6989a09327dSAlp Dener   PetscCall(TaoALMMEvaluateIterate_Private(tao));
699661095bbSAlp Dener   if (tao->eq_constrained) {
700661095bbSAlp Dener     /* compute scalar contributions */
7019566063dSJacob Faibussowitsch     PetscCall(VecDot(auglag->Ye, auglag->Ce, &yeTce));
7029566063dSJacob Faibussowitsch     PetscCall(VecDot(auglag->Ce, auglag->Ce, &ceTce));
703661095bbSAlp Dener     /* dL/dX += ye^T Ae */
7049566063dSJacob Faibussowitsch     PetscCall(MatMultTransposeAdd(auglag->Ae, auglag->Ye, auglag->LgradX, auglag->LgradX));
705a9b58be2SPaul T. Kühner     /* dL/dX += mu * ce^T Ae */
7069566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(auglag->Ae, auglag->Ce, auglag->Xwork));
707a9b58be2SPaul T. Kühner     PetscCall(VecAXPY(auglag->LgradX, auglag->mu, auglag->Xwork));
708661095bbSAlp Dener   }
709661095bbSAlp Dener   if (tao->ineq_constrained) {
710661095bbSAlp Dener     /* compute scalar contributions */
7119566063dSJacob Faibussowitsch     PetscCall(VecDot(auglag->Yi, auglag->Ci, &yiTcims));
7129566063dSJacob Faibussowitsch     PetscCall(VecDot(auglag->Ci, auglag->Ci, &cimsTcims));
713661095bbSAlp Dener     /* dL/dX += yi^T Ai */
7149566063dSJacob Faibussowitsch     PetscCall(MatMultTransposeAdd(auglag->Ai, auglag->Yi, auglag->LgradX, auglag->LgradX));
715a9b58be2SPaul T. Kühner     /* dL/dX += mu * (ci - s)^T Ai */
7169566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(auglag->Ai, auglag->Ci, auglag->Xwork));
717a9b58be2SPaul T. Kühner     PetscCall(VecAXPY(auglag->LgradX, auglag->mu, auglag->Xwork));
718661095bbSAlp Dener     /* dL/dS = -[yi + mu*(ci - s)] */
7199566063dSJacob Faibussowitsch     PetscCall(VecWAXPY(auglag->LgradS, auglag->mu, auglag->Ci, auglag->Yi));
7209566063dSJacob Faibussowitsch     PetscCall(VecScale(auglag->LgradS, -1.0));
721661095bbSAlp Dener   }
722661095bbSAlp Dener   /* combine gradient together */
7239566063dSJacob Faibussowitsch   PetscCall(TaoALMMCombinePrimal_Private(tao, auglag->LgradX, auglag->LgradS, auglag->G));
724661095bbSAlp Dener   /* compute L = f + ye^T ce + yi^T (ci - s) + 0.5*mu*||ce||^2 + 0.5*mu*||ci - s||^2 */
725661095bbSAlp Dener   auglag->Lval = auglag->fval + yeTce + yiTcims + 0.5 * auglag->mu * (ceTce + cimsTcims);
7263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
727661095bbSAlp Dener }
728661095bbSAlp Dener 
729d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoALMMSubsolverObjective_Private(Tao tao, Vec P, PetscReal *Lval, void *ctx)
730d71ae5a4SJacob Faibussowitsch {
7317721a36fSStefano Zampini   TAO_ALMM *auglag = (TAO_ALMM *)ctx;
7327721a36fSStefano Zampini 
7337721a36fSStefano Zampini   PetscFunctionBegin;
7349566063dSJacob Faibussowitsch   PetscCall(VecCopy(P, auglag->P));
7359566063dSJacob Faibussowitsch   PetscCall((*auglag->sub_obj)(auglag->parent));
7367721a36fSStefano Zampini   *Lval = auglag->Lval;
7373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7387721a36fSStefano Zampini }
7397721a36fSStefano Zampini 
740d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoALMMSubsolverObjectiveAndGradient_Private(Tao tao, Vec P, PetscReal *Lval, Vec G, void *ctx)
741d71ae5a4SJacob Faibussowitsch {
742661095bbSAlp Dener   TAO_ALMM *auglag = (TAO_ALMM *)ctx;
743661095bbSAlp Dener 
744661095bbSAlp Dener   PetscFunctionBegin;
7459566063dSJacob Faibussowitsch   PetscCall(VecCopy(P, auglag->P));
7469566063dSJacob Faibussowitsch   PetscCall((*auglag->sub_obj)(auglag->parent));
7479566063dSJacob Faibussowitsch   PetscCall(VecCopy(auglag->G, G));
748661095bbSAlp Dener   *Lval = auglag->Lval;
7493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
750661095bbSAlp Dener }
751