xref: /petsc/src/tao/constrained/impls/almm/almm.c (revision 2e956fe4fc852fabc23b437482e1fb7b77fddb0d)
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 
13661095bbSAlp Dener static PetscErrorCode TaoSolve_ALMM(Tao tao)
14661095bbSAlp Dener {
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));
259566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(auglag->Yi));
26661095bbSAlp Dener     }
27661095bbSAlp Dener     if (tao->eq_constrained) {
289566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(auglag->Ye));
29661095bbSAlp Dener     }
30661095bbSAlp Dener   }
31661095bbSAlp Dener 
32661095bbSAlp Dener   /* compute initial nonlinear Lagrangian and its derivatives */
339566063dSJacob Faibussowitsch   PetscCall((*auglag->sub_obj)(tao));
349566063dSJacob Faibussowitsch   PetscCall(TaoALMMComputeOptimalityNorms_Private(tao));
35661095bbSAlp Dener   /* print initial step and check convergence */
369566063dSJacob Faibussowitsch   PetscCall(PetscInfo(tao,"Solving with %s formulation\n",TaoALMMTypes[auglag->type]));
379566063dSJacob Faibussowitsch   PetscCall(TaoLogConvergenceHistory(tao, auglag->Lval, auglag->gnorm, auglag->cnorm, tao->ksp_its));
389566063dSJacob Faibussowitsch   PetscCall(TaoMonitor(tao, tao->niter, auglag->fval, auglag->gnorm, auglag->cnorm, 0.0));
399566063dSJacob Faibussowitsch   PetscCall((*tao->ops->convergencetest)(tao, tao->cnvP));
40661095bbSAlp Dener   /* set initial penalty factor and inner solver tolerance */
41661095bbSAlp Dener   switch (auglag->type) {
427721a36fSStefano Zampini     case TAO_ALMM_CLASSIC:
43661095bbSAlp Dener       auglag->mu = auglag->mu0;
44661095bbSAlp Dener       break;
457721a36fSStefano Zampini     case TAO_ALMM_PHR:
46661095bbSAlp Dener       auglag->cenorm = 0.0;
47661095bbSAlp Dener       if (tao->eq_constrained) {
489566063dSJacob Faibussowitsch         PetscCall(VecDot(auglag->Ce, auglag->Ce, &auglag->cenorm));
49661095bbSAlp Dener       }
50661095bbSAlp Dener       auglag->cinorm = 0.0;
51661095bbSAlp Dener       if (tao->ineq_constrained) {
529566063dSJacob Faibussowitsch         PetscCall(VecCopy(auglag->Ci, auglag->Ciwork));
539566063dSJacob Faibussowitsch         PetscCall(VecScale(auglag->Ciwork, -1.0));
549566063dSJacob Faibussowitsch         PetscCall(VecPointwiseMax(auglag->Ciwork, auglag->Cizero, auglag->Ciwork));
559566063dSJacob Faibussowitsch         PetscCall(VecDot(auglag->Ciwork, auglag->Ciwork, &auglag->cinorm));
56661095bbSAlp Dener       }
57661095bbSAlp Dener       /* determine initial penalty factor based on the balance of constraint violation and objective function value */
58661095bbSAlp Dener       auglag->mu = PetscMax(1.e-6, PetscMin(10.0, 2.0*PetscAbsReal(auglag->fval)/(auglag->cenorm + auglag->cinorm)));
59661095bbSAlp Dener       break;
60661095bbSAlp Dener     default:
61661095bbSAlp Dener       break;
62661095bbSAlp Dener   }
63661095bbSAlp Dener   auglag->gtol = auglag->gtol0;
6463a3b9bcSJacob Faibussowitsch   PetscCall(PetscInfo(tao,"Initial penalty: %.2g\n",(double)auglag->mu));
65661095bbSAlp Dener 
66661095bbSAlp Dener   /* start aug-lag outer loop */
67661095bbSAlp Dener   while (tao->reason == TAO_CONTINUE_ITERATING) {
68661095bbSAlp Dener     ++tao->niter;
69661095bbSAlp Dener     /* update subsolver tolerance */
7063a3b9bcSJacob Faibussowitsch     PetscCall(PetscInfo(tao,"Subsolver tolerance: ||G|| <= %e\n",(double)auglag->gtol));
719566063dSJacob Faibussowitsch     PetscCall(TaoSetTolerances(auglag->subsolver, auglag->gtol, 0.0, 0.0));
72661095bbSAlp Dener     /* solve the bound-constrained or unconstrained subproblem */
739566063dSJacob Faibussowitsch     PetscCall(TaoSolve(auglag->subsolver));
749566063dSJacob Faibussowitsch     PetscCall(TaoGetConvergedReason(auglag->subsolver, &reason));
75661095bbSAlp Dener     tao->ksp_its += auglag->subsolver->ksp_its;
76661095bbSAlp Dener     if (reason != TAO_CONVERGED_GATOL) {
779566063dSJacob Faibussowitsch       PetscCall(PetscInfo(tao,"Subsolver failed to converge, reason: %s\n",TaoConvergedReasons[reason]));
78661095bbSAlp Dener     }
79661095bbSAlp Dener     /* evaluate solution and test convergence */
809566063dSJacob Faibussowitsch     PetscCall((*auglag->sub_obj)(tao));
819566063dSJacob Faibussowitsch     PetscCall(TaoALMMComputeOptimalityNorms_Private(tao));
82661095bbSAlp Dener     /* decide whether to update multipliers or not */
83661095bbSAlp Dener     updated = 0.0;
84661095bbSAlp Dener     if (auglag->cnorm <= auglag->ytol) {
8563a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(tao,"Multipliers updated: ||C|| <= %e\n",(double)auglag->ytol));
86661095bbSAlp Dener       /* constraints are good, update multipliers and convergence tolerances */
87661095bbSAlp Dener       if (tao->eq_constrained) {
889566063dSJacob Faibussowitsch         PetscCall(VecAXPY(auglag->Ye, auglag->mu, auglag->Ce));
899566063dSJacob Faibussowitsch         PetscCall(VecSet(auglag->Cework, auglag->ye_max));
909566063dSJacob Faibussowitsch         PetscCall(VecPointwiseMin(auglag->Ye, auglag->Cework, auglag->Ye));
919566063dSJacob Faibussowitsch         PetscCall(VecSet(auglag->Cework, auglag->ye_min));
929566063dSJacob Faibussowitsch         PetscCall(VecPointwiseMax(auglag->Ye, auglag->Cework, auglag->Ye));
93661095bbSAlp Dener       }
94661095bbSAlp Dener       if (tao->ineq_constrained) {
959566063dSJacob Faibussowitsch         PetscCall(VecAXPY(auglag->Yi, auglag->mu, auglag->Ci));
969566063dSJacob Faibussowitsch         PetscCall(VecSet(auglag->Ciwork, auglag->yi_max));
979566063dSJacob Faibussowitsch         PetscCall(VecPointwiseMin(auglag->Yi, auglag->Ciwork, auglag->Yi));
989566063dSJacob Faibussowitsch         PetscCall(VecSet(auglag->Ciwork, auglag->yi_min));
999566063dSJacob Faibussowitsch         PetscCall(VecPointwiseMax(auglag->Yi, auglag->Ciwork, auglag->Yi));
100661095bbSAlp Dener       }
101661095bbSAlp Dener       /* tolerances are updated only for non-PHR methods */
102661095bbSAlp Dener       if (auglag->type != TAO_ALMM_PHR) {
103661095bbSAlp Dener         auglag->ytol = PetscMax(tao->catol, auglag->ytol/PetscPowReal(auglag->mu, auglag->mu_pow_good));
104661095bbSAlp Dener         auglag->gtol = PetscMax(tao->gatol, auglag->gtol/auglag->mu);
105661095bbSAlp Dener       }
106661095bbSAlp Dener       updated = 1.0;
107661095bbSAlp Dener     } else {
108661095bbSAlp Dener       /* constraints are bad, update penalty factor */
109661095bbSAlp Dener       auglag->mu = PetscMin(auglag->mu_max, auglag->mu_fac*auglag->mu);
110661095bbSAlp Dener       /* tolerances are reset only for non-PHR methods */
111661095bbSAlp Dener       if (auglag->type != TAO_ALMM_PHR) {
112661095bbSAlp Dener         auglag->ytol = PetscMax(tao->catol, 0.1/PetscPowReal(auglag->mu, auglag->mu_pow_bad));
113661095bbSAlp Dener         auglag->gtol = PetscMax(tao->gatol, 1.0/auglag->mu);
114661095bbSAlp Dener       }
11563a3b9bcSJacob Faibussowitsch       PetscCall(PetscInfo(tao,"Penalty increased: mu = %.2g\n",(double)auglag->mu));
116661095bbSAlp Dener     }
1179566063dSJacob Faibussowitsch     PetscCall(TaoLogConvergenceHistory(tao, auglag->fval, auglag->gnorm, auglag->cnorm, tao->ksp_its));
1189566063dSJacob Faibussowitsch     PetscCall(TaoMonitor(tao, tao->niter, auglag->fval, auglag->gnorm, auglag->cnorm, updated));
1199566063dSJacob Faibussowitsch     PetscCall((*tao->ops->convergencetest)(tao, tao->cnvP));
120661095bbSAlp Dener   }
121661095bbSAlp Dener 
122661095bbSAlp Dener   PetscFunctionReturn(0);
123661095bbSAlp Dener }
124661095bbSAlp Dener 
125661095bbSAlp Dener static PetscErrorCode TaoView_ALMM(Tao tao,PetscViewer viewer)
126661095bbSAlp Dener {
127661095bbSAlp Dener   TAO_ALMM       *auglag = (TAO_ALMM*)tao->data;
128661095bbSAlp Dener   PetscBool      isascii;
129661095bbSAlp Dener 
130661095bbSAlp Dener   PetscFunctionBegin;
1319566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPushTab(viewer));
1329566063dSJacob Faibussowitsch   PetscCall(TaoView(auglag->subsolver,viewer));
1339566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPopTab(viewer));
1349566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
135661095bbSAlp Dener   if (isascii) {
1369566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
1379566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "ALMM Formulation Type: %s\n", TaoALMMTypes[auglag->type]));
1389566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
139661095bbSAlp Dener   }
140661095bbSAlp Dener   PetscFunctionReturn(0);
141661095bbSAlp Dener }
142661095bbSAlp Dener 
143661095bbSAlp Dener static PetscErrorCode TaoSetUp_ALMM(Tao tao)
144661095bbSAlp Dener {
145661095bbSAlp Dener   TAO_ALMM       *auglag = (TAO_ALMM*)tao->data;
146661095bbSAlp Dener   VecType        vec_type;
147661095bbSAlp Dener   Vec            SL, SU;
148661095bbSAlp Dener   PetscBool      is_cg = PETSC_FALSE, is_lmvm = PETSC_FALSE;
149661095bbSAlp Dener 
150661095bbSAlp Dener   PetscFunctionBegin;
1513c859ba3SBarry 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.");
1523c859ba3SBarry 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.");
1539566063dSJacob Faibussowitsch   PetscCall(TaoComputeVariableBounds(tao));
154661095bbSAlp Dener   /* alias base vectors and create extras */
1559566063dSJacob Faibussowitsch   PetscCall(VecGetType(tao->solution, &vec_type));
156661095bbSAlp Dener   auglag->Px = tao->solution;
157661095bbSAlp Dener   if (!tao->gradient) { /* base gradient */
1589566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(tao->solution, &tao->gradient));
159661095bbSAlp Dener   }
160661095bbSAlp Dener   auglag->LgradX = tao->gradient;
161661095bbSAlp Dener   if (!auglag->Xwork) { /* opt var work vector */
1629566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(tao->solution, &auglag->Xwork));
163661095bbSAlp Dener   }
164661095bbSAlp Dener   if (tao->eq_constrained) {
165661095bbSAlp Dener     auglag->Ce = tao->constraints_equality;
166661095bbSAlp Dener     auglag->Ae = tao->jacobian_equality;
167661095bbSAlp Dener     if (!auglag->Ye) { /* equality multipliers */
1689566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(auglag->Ce, &auglag->Ye));
169661095bbSAlp Dener     }
170661095bbSAlp Dener     if (!auglag->Cework) {
1719566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(auglag->Ce, &auglag->Cework));
172661095bbSAlp Dener     }
173661095bbSAlp Dener   }
174661095bbSAlp Dener   if (tao->ineq_constrained) {
175661095bbSAlp Dener     auglag->Ci = tao->constraints_inequality;
176661095bbSAlp Dener     auglag->Ai = tao->jacobian_inequality;
177661095bbSAlp Dener     if (!auglag->Yi) { /* inequality multipliers */
1789566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(auglag->Ci, &auglag->Yi));
179661095bbSAlp Dener     }
180661095bbSAlp Dener     if (!auglag->Ciwork) {
1819566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(auglag->Ci, &auglag->Ciwork));
182661095bbSAlp Dener     }
183661095bbSAlp Dener     if (!auglag->Cizero) {
1849566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(auglag->Ci, &auglag->Cizero));
1859566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(auglag->Cizero));
186661095bbSAlp Dener     }
187661095bbSAlp Dener     if (!auglag->Ps) { /* slack vars */
1889566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(auglag->Ci, &auglag->Ps));
189661095bbSAlp Dener     }
190661095bbSAlp Dener     if (!auglag->LgradS) { /* slack component of Lagrangian gradient */
1919566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(auglag->Ci, &auglag->LgradS));
192661095bbSAlp Dener     }
193661095bbSAlp Dener     /* create vector for combined primal space and the associated communication objects */
194661095bbSAlp Dener     if (!auglag->P) {
1959566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(2, &auglag->Parr));
196661095bbSAlp Dener       auglag->Parr[0] = auglag->Px; auglag->Parr[1] = auglag->Ps;
1979566063dSJacob Faibussowitsch       PetscCall(VecConcatenate(2, auglag->Parr, &auglag->P, &auglag->Pis));
1989566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(2, &auglag->Pscatter));
1999566063dSJacob Faibussowitsch       PetscCall(VecScatterCreate(auglag->P, auglag->Pis[0], auglag->Px, NULL, &auglag->Pscatter[0]));
2009566063dSJacob Faibussowitsch       PetscCall(VecScatterCreate(auglag->P, auglag->Pis[1], auglag->Ps, NULL, &auglag->Pscatter[1]));
201661095bbSAlp Dener     }
202661095bbSAlp Dener     if (tao->eq_constrained) {
203661095bbSAlp Dener       /* create vector for combined dual space and the associated communication objects */
204661095bbSAlp Dener       if (!auglag->Y) {
2059566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(2, &auglag->Yarr));
206661095bbSAlp Dener         auglag->Yarr[0] = auglag->Ye; auglag->Yarr[1] = auglag->Yi;
2079566063dSJacob Faibussowitsch         PetscCall(VecConcatenate(2, auglag->Yarr, &auglag->Y, &auglag->Yis));
2089566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(2, &auglag->Yscatter));
2099566063dSJacob Faibussowitsch         PetscCall(VecScatterCreate(auglag->Y, auglag->Yis[0], auglag->Ye, NULL, &auglag->Yscatter[0]));
2109566063dSJacob Faibussowitsch         PetscCall(VecScatterCreate(auglag->Y, auglag->Yis[1], auglag->Yi, NULL, &auglag->Yscatter[1]));
211661095bbSAlp Dener       }
212661095bbSAlp Dener       if (!auglag->C) {
2139566063dSJacob Faibussowitsch         PetscCall(VecDuplicate(auglag->Y, &auglag->C));
214661095bbSAlp Dener       }
215661095bbSAlp Dener     } else {
216661095bbSAlp Dener       if (!auglag->C) {
217661095bbSAlp Dener         auglag->C = auglag->Ci;
218661095bbSAlp Dener       }
219661095bbSAlp Dener       if (!auglag->Y) {
220661095bbSAlp Dener         auglag->Y = auglag->Yi;
221661095bbSAlp Dener       }
222661095bbSAlp Dener     }
223661095bbSAlp Dener   } else {
224661095bbSAlp Dener     if (!auglag->P) {
225661095bbSAlp Dener       auglag->P = auglag->Px;
226661095bbSAlp Dener     }
227661095bbSAlp Dener     if (!auglag->G) {
228661095bbSAlp Dener       auglag->G = auglag->LgradX;
229661095bbSAlp Dener     }
230661095bbSAlp Dener     if (!auglag->C) {
231661095bbSAlp Dener       auglag->C = auglag->Ce;
232661095bbSAlp Dener     }
233661095bbSAlp Dener     if (!auglag->Y) {
234661095bbSAlp Dener       auglag->Y = auglag->Ye;
235661095bbSAlp Dener     }
236661095bbSAlp Dener   }
237661095bbSAlp Dener   /* initialize parameters */
238661095bbSAlp Dener   if (auglag->type == TAO_ALMM_PHR) {
239661095bbSAlp Dener     auglag->mu_fac = 10.0;
240661095bbSAlp Dener     auglag->yi_min = 0.0;
241661095bbSAlp Dener     auglag->ytol0 = 0.5;
242661095bbSAlp Dener     auglag->gtol0 = tao->gatol;
243661095bbSAlp Dener     if (tao->gatol_changed && tao->catol_changed) {
2449566063dSJacob Faibussowitsch       PetscCall(PetscInfo(tao,"TAOALMM with PHR: different gradient and constraint tolerances are not supported, setting catol = gatol\n"));
245661095bbSAlp Dener       tao->catol = tao->gatol;
246661095bbSAlp Dener     }
247661095bbSAlp Dener   }
248661095bbSAlp Dener   /* set the Lagrangian formulation type for the subsolver */
249661095bbSAlp Dener   switch (auglag->type) {
2507721a36fSStefano Zampini     case TAO_ALMM_CLASSIC:
251661095bbSAlp Dener       auglag->sub_obj = TaoALMMComputeAugLagAndGradient_Private;
252661095bbSAlp Dener       break;
2537721a36fSStefano Zampini     case TAO_ALMM_PHR:
254661095bbSAlp Dener       auglag->sub_obj = TaoALMMComputePHRLagAndGradient_Private;
255661095bbSAlp Dener       break;
256661095bbSAlp Dener     default:
257661095bbSAlp Dener       break;
258661095bbSAlp Dener   }
259661095bbSAlp Dener   /* set up the subsolver */
2609566063dSJacob Faibussowitsch   PetscCall(TaoSetSolution(auglag->subsolver, auglag->P));
2619566063dSJacob Faibussowitsch   PetscCall(TaoSetObjective(auglag->subsolver, TaoALMMSubsolverObjective_Private, (void*)auglag));
2629566063dSJacob Faibussowitsch   PetscCall(TaoSetObjectiveAndGradient(auglag->subsolver, NULL, TaoALMMSubsolverObjectiveAndGradient_Private, (void*)auglag));
263661095bbSAlp Dener   if (tao->bounded) {
264661095bbSAlp Dener     /* make sure that the subsolver is a bound-constrained method */
2659566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)auglag->subsolver, TAOCG, &is_cg));
2669566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)auglag->subsolver, TAOLMVM, &is_lmvm));
267661095bbSAlp Dener     if (is_cg) {
2689566063dSJacob Faibussowitsch       PetscCall(TaoSetType(auglag->subsolver, TAOBNCG));
2699566063dSJacob Faibussowitsch       PetscCall(PetscInfo(tao,"TAOCG detected for bound-constrained problem, switching to TAOBNCG instead."));
270661095bbSAlp Dener     }
271661095bbSAlp Dener     if (is_lmvm) {
2729566063dSJacob Faibussowitsch       PetscCall(TaoSetType(auglag->subsolver, TAOBQNLS));
2739566063dSJacob Faibussowitsch       PetscCall(PetscInfo(tao,"TAOLMVM detected for bound-constrained problem, switching to TAOBQNLS instead."));
274661095bbSAlp Dener     }
275661095bbSAlp Dener     /* create lower and upper bound clone vectors for subsolver */
276661095bbSAlp Dener     if (!auglag->PL) {
2779566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(auglag->P, &auglag->PL));
278661095bbSAlp Dener     }
279661095bbSAlp Dener     if (!auglag->PU) {
2809566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(auglag->P, &auglag->PU));
281661095bbSAlp Dener     }
282661095bbSAlp Dener     if (tao->ineq_constrained) {
283661095bbSAlp Dener       /* create lower and upper bounds for slack, set lower to 0 */
2849566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(auglag->Ci, &SL));
2859566063dSJacob Faibussowitsch       PetscCall(VecSet(SL, 0.0));
2869566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(auglag->Ci, &SU));
2879566063dSJacob Faibussowitsch       PetscCall(VecSet(SU, PETSC_INFINITY));
288661095bbSAlp Dener       /* combine opt var bounds with slack bounds */
2899566063dSJacob Faibussowitsch       PetscCall(TaoALMMCombinePrimal_Private(tao, tao->XL, SL, auglag->PL));
2909566063dSJacob Faibussowitsch       PetscCall(TaoALMMCombinePrimal_Private(tao, tao->XU, SU, auglag->PU));
291661095bbSAlp Dener       /* destroy work vectors */
2929566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&SL));
2939566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&SU));
294661095bbSAlp Dener     } else {
295661095bbSAlp Dener       /* no inequality constraints, just copy bounds into the subsolver */
2969566063dSJacob Faibussowitsch       PetscCall(VecCopy(tao->XL, auglag->PL));
2979566063dSJacob Faibussowitsch       PetscCall(VecCopy(tao->XU, auglag->PU));
298661095bbSAlp Dener     }
2999566063dSJacob Faibussowitsch     PetscCall(TaoSetVariableBounds(auglag->subsolver, auglag->PL, auglag->PU));
300661095bbSAlp Dener   }
3019566063dSJacob Faibussowitsch   PetscCall(TaoSetUp(auglag->subsolver));
302661095bbSAlp Dener   auglag->G = auglag->subsolver->gradient;
303661095bbSAlp Dener 
304661095bbSAlp Dener   PetscFunctionReturn(0);
305661095bbSAlp Dener }
306661095bbSAlp Dener 
307661095bbSAlp Dener static PetscErrorCode TaoDestroy_ALMM(Tao tao)
308661095bbSAlp Dener {
309661095bbSAlp Dener   TAO_ALMM       *auglag = (TAO_ALMM*)tao->data;
310661095bbSAlp Dener 
311661095bbSAlp Dener   PetscFunctionBegin;
3129566063dSJacob Faibussowitsch   PetscCall(TaoDestroy(&auglag->subsolver));
313661095bbSAlp Dener   if (tao->setupcalled) {
3149566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&auglag->Xwork));              /* opt work */
315661095bbSAlp Dener     if (tao->eq_constrained) {
3169566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&auglag->Ye));               /* equality multipliers */
3179566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&auglag->Cework));           /* equality work vector */
318661095bbSAlp Dener     }
319661095bbSAlp Dener     if (tao->ineq_constrained) {
3209566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&auglag->Ps));               /* slack vars */
321661095bbSAlp Dener       auglag->Parr[0] = NULL;                                     /* clear pointer to tao->solution, will be destroyed by TaoDestroy() shell */
3229566063dSJacob Faibussowitsch       PetscCall(PetscFree(auglag->Parr));               /* array of primal vectors */
3239566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&auglag->LgradS));           /* slack grad */
3249566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&auglag->Cizero));           /* zero vector for pointwise max */
3259566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&auglag->Yi));               /* inequality multipliers */
3269566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&auglag->Ciwork));           /* inequality work vector */
3279566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&auglag->P));                /* combo primal */
3289566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&auglag->Pis[0]));            /* index set for X inside P */
3299566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&auglag->Pis[1]));            /* index set for S inside P */
3309566063dSJacob Faibussowitsch       PetscCall(PetscFree(auglag->Pis));                /* array of P index sets */
3319566063dSJacob Faibussowitsch       PetscCall(VecScatterDestroy(&auglag->Pscatter[0]));
3329566063dSJacob Faibussowitsch       PetscCall(VecScatterDestroy(&auglag->Pscatter[1]));
3339566063dSJacob Faibussowitsch       PetscCall(PetscFree(auglag->Pscatter));
334661095bbSAlp Dener       if (tao->eq_constrained) {
3359566063dSJacob Faibussowitsch         PetscCall(VecDestroy(&auglag->Y));              /* combo multipliers */
3369566063dSJacob Faibussowitsch         PetscCall(PetscFree(auglag->Yarr));             /* array of dual vectors */
3379566063dSJacob Faibussowitsch         PetscCall(VecDestroy(&auglag->C));              /* combo constraints */
3389566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&auglag->Yis[0]));          /* index set for Ye inside Y */
3399566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&auglag->Yis[1]));          /* index set for Yi inside Y */
3409566063dSJacob Faibussowitsch         PetscCall(PetscFree(auglag->Yis));
3419566063dSJacob Faibussowitsch         PetscCall(VecScatterDestroy(&auglag->Yscatter[0]));
3429566063dSJacob Faibussowitsch         PetscCall(VecScatterDestroy(&auglag->Yscatter[1]));
3439566063dSJacob Faibussowitsch         PetscCall(PetscFree(auglag->Yscatter));
344661095bbSAlp Dener       }
345661095bbSAlp Dener     }
346661095bbSAlp Dener     if (tao->bounded) {
3479566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&auglag->PL));                /* lower bounds for subsolver */
3489566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&auglag->PU));                /* upper bounds for subsolver */
349661095bbSAlp Dener     }
350661095bbSAlp Dener   }
351*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoALMMGetType_C",NULL));
352*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoALMMSetType_C",NULL));
353*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoALMMGetSubsolver_C",NULL));
354*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoALMMSetSubsolver_C",NULL));
355*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoALMMGetMultipliers_C",NULL));
356*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoALMMSetMultipliers_C",NULL));
357*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoALMMGetPrimalIS_C",NULL));
358*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoALMMGetDualIS_C",NULL));
3599566063dSJacob Faibussowitsch   PetscCall(PetscFree(tao->data));
360661095bbSAlp Dener   PetscFunctionReturn(0);
361661095bbSAlp Dener }
362661095bbSAlp Dener 
363661095bbSAlp Dener static PetscErrorCode TaoSetFromOptions_ALMM(PetscOptionItems *PetscOptionsObject,Tao tao)
364661095bbSAlp Dener {
365661095bbSAlp Dener   TAO_ALMM       *auglag = (TAO_ALMM*)tao->data;
366661095bbSAlp Dener   PetscInt       i;
367661095bbSAlp Dener 
368661095bbSAlp Dener   PetscFunctionBegin;
36971075aafSSatish Balay   PetscOptionsHeadBegin(PetscOptionsObject,"Augmented Lagrangian multiplier method solves problems with general constraints by converting them into a sequence of unconstrained problems.");
3709566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_almm_mu_init","initial penalty parameter","",auglag->mu0,&auglag->mu0,NULL));
3719566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_almm_mu_factor","increase factor for the penalty parameter","",auglag->mu_fac,&auglag->mu_fac,NULL));
3729566063dSJacob 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));
3739566063dSJacob 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));
3749566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_almm_mu_max","maximum safeguard for penalty parameter updates","",auglag->mu_max,&auglag->mu_max,NULL));
3759566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_almm_ye_min","minimum safeguard for equality multiplier updates","",auglag->ye_min,&auglag->ye_min,NULL));
3769566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_almm_ye_max","maximum safeguard for equality multipliers updates","",auglag->ye_max,&auglag->ye_max,NULL));
3779566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_almm_yi_min","minimum safeguard for inequality multipliers updates","",auglag->yi_min,&auglag->yi_min,NULL));
3789566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tao_almm_yi_max","maximum safeguard for inequality multipliers updates","",auglag->yi_max,&auglag->yi_max,NULL));
3799566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-tao_almm_type","augmented Lagrangian formulation type for the subproblem","TaoALMMType",TaoALMMTypes,(PetscEnum)auglag->type,(PetscEnum*)&auglag->type,NULL));
380d0609cedSBarry Smith   PetscOptionsHeadEnd();
3819566063dSJacob Faibussowitsch   PetscCall(TaoSetOptionsPrefix(auglag->subsolver,((PetscObject)tao)->prefix));
3829566063dSJacob Faibussowitsch   PetscCall(TaoAppendOptionsPrefix(auglag->subsolver,"tao_almm_subsolver_"));
3839566063dSJacob Faibussowitsch   PetscCall(TaoSetFromOptions(auglag->subsolver));
384661095bbSAlp Dener   for (i=0; i<tao->numbermonitors; i++) {
3859566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)tao->monitorcontext[i]));
3869566063dSJacob Faibussowitsch     PetscCall(TaoSetMonitor(auglag->subsolver, tao->monitor[i], tao->monitorcontext[i], tao->monitordestroy[i]));
387661095bbSAlp Dener     if (tao->monitor[i] == TaoMonitorDefault || tao->monitor[i] == TaoDefaultCMonitor || tao->monitor[i] == TaoDefaultGMonitor || tao->monitor[i] == TaoDefaultSMonitor) {
388661095bbSAlp Dener       auglag->info = PETSC_TRUE;
389661095bbSAlp Dener     }
390661095bbSAlp Dener   }
391661095bbSAlp Dener   PetscFunctionReturn(0);
392661095bbSAlp Dener }
393661095bbSAlp Dener 
394661095bbSAlp Dener /* -------------------------------------------------------- */
395661095bbSAlp Dener 
396661095bbSAlp Dener /*MC
397661095bbSAlp Dener   TaoALMM - Augmented Lagrangian multiplier method for solving nonlinear optimization problems with general constraints.
398661095bbSAlp Dener 
399661095bbSAlp Dener   Options Database Keys:
400661095bbSAlp Dener + -tao_almm_mu_init <real>       - initial penalty parameter (default: 10.)
401661095bbSAlp Dener . -tao_almm_mu_factor <real>     - increase factor for the penalty parameter (default: 100.)
402661095bbSAlp Dener . -tao_almm_mu_max <real>        - maximum safeguard for penalty parameter updates (default: 1.e20)
403661095bbSAlp Dener . -tao_almm_mu_power_good <real> - exponential for penalty parameter when multiplier update is accepted (default: 0.9)
404661095bbSAlp Dener . -tao_almm_mu_power_bad <real>  - exponential for penalty parameter when multiplier update is rejected (default: 0.1)
405661095bbSAlp Dener . -tao_almm_ye_min <real>        - minimum safeguard for equality multiplier updates (default: -1.e20)
406661095bbSAlp Dener . -tao_almm_ye_max <real>        - maximum safeguard for equality multiplier updates (default: 1.e20)
407661095bbSAlp Dener . -tao_almm_yi_min <real>        - minimum safeguard for inequality multiplier updates (default: -1.e20)
408661095bbSAlp Dener . -tao_almm_yi_max <real>        - maximum safeguard for inequality multiplier updates (default: 1.e20)
409661095bbSAlp Dener - -tao_almm_type <classic,phr>   - change formulation of the augmented Lagrangian merit function for the subproblem (default: classic)
410661095bbSAlp Dener 
411661095bbSAlp Dener   Level: beginner
412661095bbSAlp Dener 
413661095bbSAlp Dener   Notes:
414661095bbSAlp Dener   This method converts a constrained problem into a sequence of unconstrained problems via the augmented
415661095bbSAlp Dener   Lagrangian merit function. Bound constraints are pushed down to the subproblem without any modifications.
416661095bbSAlp Dener 
417661095bbSAlp Dener   Two formulations are offered for the subproblem: canonical Hestenes-Powell augmented Lagrangian with slack
418661095bbSAlp Dener   variables for inequality constraints, and a slack-less Powell-Hestenes-Rockafellar (PHR) formulation utilizing a
419661095bbSAlp Dener   pointwise max() penalty on inequality constraints. The canonical augmented Lagrangian formulation typically
420661095bbSAlp Dener   converges faster for most problems. However, PHR may be desirable for problems featuring a large number
421661095bbSAlp Dener   of inequality constraints because it avoids inflating the size of the subproblem with slack variables.
422661095bbSAlp Dener 
423661095bbSAlp Dener   The subproblem is solved using a nested first-order TAO solver. The user can retrieve a pointer to
424661095bbSAlp Dener   the subsolver via TaoALMMGetSubsolver() or pass command line arguments to it using the
425661095bbSAlp Dener   "-tao_almm_subsolver_" prefix. Currently, TaoALMM does not support second-order methods for the
426661095bbSAlp Dener   subproblem. It is also highly recommended that the subsolver chosen by the user utilize a trust-region
427661095bbSAlp Dener   strategy for globalization (default: TAOBQNKTR) especially if the outer problem features bound constraints.
428661095bbSAlp Dener 
429661095bbSAlp Dener .vb
430661095bbSAlp Dener   while unconverged
431661095bbSAlp Dener     solve argmin_x L(x) s.t. l <= x <= u
432661095bbSAlp Dener     if ||c|| <= y_tol
433661095bbSAlp Dener       if ||c|| <= c_tol && ||Lgrad|| <= g_tol:
434661095bbSAlp Dener         problem converged, return solution
435661095bbSAlp Dener       else
436661095bbSAlp Dener         constraints sufficiently improved
437661095bbSAlp Dener         update multipliers and tighten tolerances
438661095bbSAlp Dener       endif
439661095bbSAlp Dener     else
440661095bbSAlp Dener       constraints did not improve
441661095bbSAlp Dener       update penalty and loosen tolerances
442661095bbSAlp Dener     endif
443661095bbSAlp Dener   endwhile
444661095bbSAlp Dener .ve
445661095bbSAlp Dener 
446db781477SPatrick Sanan .seealso: `TaoALMMGetType()`, `TaoALMMSetType()`, `TaoALMMSetSubsolver()`, `TaoALMMGetSubsolver()`,
447db781477SPatrick Sanan           `TaoALMMGetMultipliers()`, `TaoALMMSetMultipliers()`, `TaoALMMGetPrimalIS()`, `TaoALMMGetDualIS()`
448661095bbSAlp Dener M*/
449661095bbSAlp Dener PETSC_EXTERN PetscErrorCode TaoCreate_ALMM(Tao tao)
450661095bbSAlp Dener {
451661095bbSAlp Dener   TAO_ALMM       *auglag;
452661095bbSAlp Dener 
453661095bbSAlp Dener   PetscFunctionBegin;
4549566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(tao, &auglag));
455661095bbSAlp Dener 
456661095bbSAlp Dener   tao->ops->destroy        = TaoDestroy_ALMM;
457661095bbSAlp Dener   tao->ops->setup          = TaoSetUp_ALMM;
458661095bbSAlp Dener   tao->ops->setfromoptions = TaoSetFromOptions_ALMM;
459661095bbSAlp Dener   tao->ops->view           = TaoView_ALMM;
460661095bbSAlp Dener   tao->ops->solve          = TaoSolve_ALMM;
461661095bbSAlp Dener 
462661095bbSAlp Dener   tao->gatol = 1.e-5;
463661095bbSAlp Dener   tao->grtol = 0.0;
464661095bbSAlp Dener   tao->gttol = 0.0;
465661095bbSAlp Dener   tao->catol = 1.e-5;
466661095bbSAlp Dener   tao->crtol = 0.0;
467661095bbSAlp Dener 
468661095bbSAlp Dener   tao->data           = (void*)auglag;
469661095bbSAlp Dener   auglag->parent      = tao;
470661095bbSAlp Dener   auglag->mu0         = 10.0;
471661095bbSAlp Dener   auglag->mu          = auglag->mu0;
472661095bbSAlp Dener   auglag->mu_fac      = 10.0;
473661095bbSAlp Dener   auglag->mu_max      = PETSC_INFINITY;
474661095bbSAlp Dener   auglag->mu_pow_good = 0.9;
475661095bbSAlp Dener   auglag->mu_pow_bad  = 0.1;
476661095bbSAlp Dener   auglag->ye_min      = PETSC_NINFINITY;
477661095bbSAlp Dener   auglag->ye_max      = PETSC_INFINITY;
478661095bbSAlp Dener   auglag->yi_min      = PETSC_NINFINITY;
479661095bbSAlp Dener   auglag->yi_max      = PETSC_INFINITY;
480661095bbSAlp Dener   auglag->ytol0       = 0.1/PetscPowReal(auglag->mu0, auglag->mu_pow_bad);
481661095bbSAlp Dener   auglag->ytol        = auglag->ytol0;
482661095bbSAlp Dener   auglag->gtol0       = 1.0/auglag->mu0;
483661095bbSAlp Dener   auglag->gtol        = auglag->gtol0;
484661095bbSAlp Dener 
485661095bbSAlp Dener   auglag->sub_obj     = TaoALMMComputeAugLagAndGradient_Private;
486661095bbSAlp Dener   auglag->type        = TAO_ALMM_CLASSIC;
487661095bbSAlp Dener   auglag->info        = PETSC_FALSE;
488661095bbSAlp Dener 
4899566063dSJacob Faibussowitsch   PetscCall(TaoCreate(PetscObjectComm((PetscObject)tao),&auglag->subsolver));
4909566063dSJacob Faibussowitsch   PetscCall(TaoSetType(auglag->subsolver, TAOBQNKTR));
4919566063dSJacob Faibussowitsch   PetscCall(TaoSetTolerances(auglag->subsolver, auglag->gtol, 0.0, 0.0));
4929566063dSJacob Faibussowitsch   PetscCall(TaoSetMaximumIterations(auglag->subsolver, 1000));
4939566063dSJacob Faibussowitsch   PetscCall(TaoSetMaximumFunctionEvaluations(auglag->subsolver, 10000));
4949566063dSJacob Faibussowitsch   PetscCall(TaoSetFunctionLowerBound(auglag->subsolver, PETSC_NINFINITY));
4959566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)auglag->subsolver,(PetscObject)tao,1));
496661095bbSAlp Dener 
4979566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetType_C", TaoALMMGetType_Private));
4989566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMSetType_C", TaoALMMSetType_Private));
4999566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetSubsolver_C", TaoALMMGetSubsolver_Private));
5009566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMSetSubsolver_C", TaoALMMSetSubsolver_Private));
5019566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetMultipliers_C", TaoALMMGetMultipliers_Private));
5029566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMSetMultipliers_C", TaoALMMSetMultipliers_Private));
5039566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetPrimalIS_C", TaoALMMGetPrimalIS_Private));
5049566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetDualIS_C", TaoALMMGetDualIS_Private));
505661095bbSAlp Dener   PetscFunctionReturn(0);
506661095bbSAlp Dener }
507661095bbSAlp Dener 
508661095bbSAlp Dener static PetscErrorCode TaoALMMCombinePrimal_Private(Tao tao, Vec X, Vec S, Vec P)
509661095bbSAlp Dener {
510661095bbSAlp Dener   TAO_ALMM     *auglag = (TAO_ALMM*)tao->data;
511661095bbSAlp Dener 
512661095bbSAlp Dener   PetscFunctionBegin;
513661095bbSAlp Dener   if (tao->ineq_constrained) {
5149566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(auglag->Pscatter[0], X, P, INSERT_VALUES, SCATTER_REVERSE));
5159566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(auglag->Pscatter[0], X, P, INSERT_VALUES, SCATTER_REVERSE));
5169566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(auglag->Pscatter[1], S, P, INSERT_VALUES, SCATTER_REVERSE));
5179566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(auglag->Pscatter[1], S, P, INSERT_VALUES, SCATTER_REVERSE));
518661095bbSAlp Dener   } else {
5199566063dSJacob Faibussowitsch     PetscCall(VecCopy(X, P));
520661095bbSAlp Dener   }
521661095bbSAlp Dener   PetscFunctionReturn(0);
522661095bbSAlp Dener }
523661095bbSAlp Dener 
524661095bbSAlp Dener static PetscErrorCode TaoALMMCombineDual_Private(Tao tao, Vec EQ, Vec IN, Vec Y)
525661095bbSAlp Dener {
526661095bbSAlp Dener   TAO_ALMM     *auglag = (TAO_ALMM*)tao->data;
527661095bbSAlp Dener 
528661095bbSAlp Dener   PetscFunctionBegin;
529661095bbSAlp Dener   if (tao->eq_constrained) {
530661095bbSAlp Dener     if (tao->ineq_constrained) {
5319566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(auglag->Yscatter[0], EQ, Y, INSERT_VALUES, SCATTER_REVERSE));
5329566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(auglag->Yscatter[0], EQ, Y, INSERT_VALUES, SCATTER_REVERSE));
5339566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(auglag->Yscatter[1], IN, Y, INSERT_VALUES, SCATTER_REVERSE));
5349566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(auglag->Yscatter[1], IN, Y, INSERT_VALUES, SCATTER_REVERSE));
535661095bbSAlp Dener     } else {
5369566063dSJacob Faibussowitsch       PetscCall(VecCopy(EQ, Y));
537661095bbSAlp Dener     }
538661095bbSAlp Dener   } else {
5399566063dSJacob Faibussowitsch     PetscCall(VecCopy(IN, Y));
540661095bbSAlp Dener   }
541661095bbSAlp Dener   PetscFunctionReturn(0);
542661095bbSAlp Dener }
543661095bbSAlp Dener 
544661095bbSAlp Dener static PetscErrorCode TaoALMMSplitPrimal_Private(Tao tao, Vec P, Vec X, Vec S)
545661095bbSAlp Dener {
546661095bbSAlp Dener   TAO_ALMM     *auglag = (TAO_ALMM*)tao->data;
547661095bbSAlp Dener 
548661095bbSAlp Dener   PetscFunctionBegin;
549661095bbSAlp Dener   if (tao->ineq_constrained) {
5509566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(auglag->Pscatter[0], P, X, INSERT_VALUES, SCATTER_FORWARD));
5519566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(auglag->Pscatter[0], P, X, INSERT_VALUES, SCATTER_FORWARD));
5529566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(auglag->Pscatter[1], P, S, INSERT_VALUES, SCATTER_FORWARD));
5539566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(auglag->Pscatter[1], P, S, INSERT_VALUES, SCATTER_FORWARD));
554661095bbSAlp Dener   } else {
5559566063dSJacob Faibussowitsch     PetscCall(VecCopy(P, X));
556661095bbSAlp Dener   }
557661095bbSAlp Dener   PetscFunctionReturn(0);
558661095bbSAlp Dener }
559661095bbSAlp Dener 
560661095bbSAlp Dener /* this assumes that the latest constraints are stored in Ce and Ci, and also combined in C */
561661095bbSAlp Dener static PetscErrorCode TaoALMMComputeOptimalityNorms_Private(Tao tao)
562661095bbSAlp Dener {
563661095bbSAlp Dener   TAO_ALMM     *auglag = (TAO_ALMM*)tao->data;
564661095bbSAlp Dener 
565661095bbSAlp Dener   PetscFunctionBegin;
566661095bbSAlp Dener   /* if bounded, project the gradient */
567661095bbSAlp Dener   if (tao->bounded) {
5689566063dSJacob Faibussowitsch     PetscCall(VecBoundGradientProjection(auglag->LgradX, auglag->Px, tao->XL, tao->XU, auglag->LgradX));
569661095bbSAlp Dener   }
570661095bbSAlp Dener   if (auglag->type == TAO_ALMM_PHR) {
5719566063dSJacob Faibussowitsch     PetscCall(VecNorm(auglag->LgradX, NORM_INFINITY, &auglag->gnorm));
572661095bbSAlp Dener     auglag->cenorm = 0.0;
573661095bbSAlp Dener     if (tao->eq_constrained) {
5749566063dSJacob Faibussowitsch       PetscCall(VecNorm(auglag->Ce, NORM_INFINITY, &auglag->cenorm));
575661095bbSAlp Dener     }
576661095bbSAlp Dener     auglag->cinorm = 0.0;
577661095bbSAlp Dener     if (tao->ineq_constrained) {
5789566063dSJacob Faibussowitsch       PetscCall(VecCopy(auglag->Yi, auglag->Ciwork));
5799566063dSJacob Faibussowitsch       PetscCall(VecScale(auglag->Ciwork, -1.0/auglag->mu));
5809566063dSJacob Faibussowitsch       PetscCall(VecPointwiseMax(auglag->Ciwork, auglag->Ci, auglag->Ciwork));
5819566063dSJacob Faibussowitsch       PetscCall(VecNorm(auglag->Ciwork, NORM_INFINITY, &auglag->cinorm));
582661095bbSAlp Dener     }
583661095bbSAlp Dener     auglag->cnorm_old = auglag->cnorm;
584661095bbSAlp Dener     auglag->cnorm = PetscMax(auglag->cenorm, auglag->cinorm);
585661095bbSAlp Dener     auglag->ytol = auglag->ytol0 * auglag->cnorm_old;
586661095bbSAlp Dener   } else {
5879566063dSJacob Faibussowitsch     PetscCall(VecNorm(auglag->LgradX, NORM_2, &auglag->gnorm));
5889566063dSJacob Faibussowitsch     PetscCall(VecNorm(auglag->C, NORM_2, &auglag->cnorm));
589661095bbSAlp Dener   }
590661095bbSAlp Dener   PetscFunctionReturn(0);
591661095bbSAlp Dener }
592661095bbSAlp Dener 
593661095bbSAlp Dener static PetscErrorCode TaoALMMEvaluateIterate_Private(Tao tao, Vec P)
594661095bbSAlp Dener {
595661095bbSAlp Dener   TAO_ALMM       *auglag = (TAO_ALMM*)tao->data;
596661095bbSAlp Dener 
597661095bbSAlp Dener   PetscFunctionBegin;
598661095bbSAlp Dener   /* split solution into primal and slack components */
5999566063dSJacob Faibussowitsch   PetscCall(TaoALMMSplitPrimal_Private(tao, auglag->P, auglag->Px, auglag->Ps));
600661095bbSAlp Dener 
601661095bbSAlp Dener   /* compute f, df/dx and the constraints */
6029566063dSJacob Faibussowitsch   PetscCall(TaoComputeObjectiveAndGradient(tao, auglag->Px, &auglag->fval, auglag->LgradX));
603661095bbSAlp Dener   if (tao->eq_constrained) {
6049566063dSJacob Faibussowitsch     PetscCall(TaoComputeEqualityConstraints(tao, auglag->Px, auglag->Ce));
6059566063dSJacob Faibussowitsch     PetscCall(TaoComputeJacobianEquality(tao, auglag->Px, auglag->Ae, auglag->Ae));
606661095bbSAlp Dener   }
607661095bbSAlp Dener   if (tao->ineq_constrained) {
6089566063dSJacob Faibussowitsch     PetscCall(TaoComputeInequalityConstraints(tao, auglag->Px, auglag->Ci));
6099566063dSJacob Faibussowitsch     PetscCall(TaoComputeJacobianInequality(tao, auglag->Px, auglag->Ai, auglag->Ai));
610661095bbSAlp Dener     switch (auglag->type) {
6117721a36fSStefano Zampini       case TAO_ALMM_CLASSIC:
612661095bbSAlp Dener         /* classic formulation converts inequality to equality constraints via slack variables */
6139566063dSJacob Faibussowitsch         PetscCall(VecAXPY(auglag->Ci, -1.0, auglag->Ps));
614661095bbSAlp Dener         break;
6157721a36fSStefano Zampini       case TAO_ALMM_PHR:
616661095bbSAlp Dener         /* PHR is based on Ci <= 0 while TAO defines Ci >= 0 so we hit it with a negative sign */
6179566063dSJacob Faibussowitsch         PetscCall(VecScale(auglag->Ci, -1.0));
6189566063dSJacob Faibussowitsch         PetscCall(MatScale(auglag->Ai, -1.0));
619661095bbSAlp Dener         break;
620661095bbSAlp Dener       default:
621661095bbSAlp Dener         break;
622661095bbSAlp Dener     }
623661095bbSAlp Dener   }
624661095bbSAlp Dener   /* combine constraints into one vector */
6259566063dSJacob Faibussowitsch   PetscCall(TaoALMMCombineDual_Private(tao, auglag->Ce, auglag->Ci, auglag->C));
626661095bbSAlp Dener   PetscFunctionReturn(0);
627661095bbSAlp Dener }
628661095bbSAlp Dener 
629661095bbSAlp Dener /*
630661095bbSAlp Dener Lphr = f + 0.5*mu*[ (Ce + Ye/mu)^T (Ce + Ye/mu) + pmin(0, Ci + Yi/mu)^T pmin(0, Ci + Yi/mu)]
631661095bbSAlp Dener 
632661095bbSAlp Dener dLphr/dX = dF/dX + mu*[ (Ce + Ye/mu)^T Ae + pmin(0, Ci + Yi/mu)^T Ai]
633661095bbSAlp Dener 
634661095bbSAlp Dener dLphr/dS = 0
635661095bbSAlp Dener */
636661095bbSAlp Dener static PetscErrorCode TaoALMMComputePHRLagAndGradient_Private(Tao tao)
637661095bbSAlp Dener {
638661095bbSAlp Dener   TAO_ALMM       *auglag = (TAO_ALMM*)tao->data;
639661095bbSAlp Dener   PetscReal      eq_norm=0.0, ineq_norm=0.0;
640661095bbSAlp Dener 
641661095bbSAlp Dener   PetscFunctionBegin;
6429566063dSJacob Faibussowitsch   PetscCall(TaoALMMEvaluateIterate_Private(tao, auglag->P));
643661095bbSAlp Dener   if (tao->eq_constrained) {
644661095bbSAlp Dener     /* Ce_work = mu*(Ce + Ye/mu) */
6459566063dSJacob Faibussowitsch     PetscCall(VecWAXPY(auglag->Cework, 1.0/auglag->mu, auglag->Ye, auglag->Ce));
6469566063dSJacob Faibussowitsch     PetscCall(VecDot(auglag->Cework, auglag->Cework, &eq_norm)); /* contribution to scalar Lagrangian */
6479566063dSJacob Faibussowitsch     PetscCall(VecScale(auglag->Cework, auglag->mu));
648661095bbSAlp Dener     /* dL/dX += mu*(Ce + Ye/mu)^T Ae */
6499566063dSJacob Faibussowitsch     PetscCall(MatMultTransposeAdd(auglag->Ae, auglag->Cework, auglag->LgradX, auglag->LgradX));
650661095bbSAlp Dener   }
651661095bbSAlp Dener   if (tao->ineq_constrained) {
652661095bbSAlp Dener     /* Ci_work = mu * pmax(0, Ci + Yi/mu) where pmax() is pointwise max() */
6539566063dSJacob Faibussowitsch     PetscCall(VecWAXPY(auglag->Ciwork, 1.0/auglag->mu, auglag->Yi, auglag->Ci));
6549566063dSJacob Faibussowitsch     PetscCall(VecPointwiseMax(auglag->Ciwork, auglag->Cizero, auglag->Ciwork));
6559566063dSJacob Faibussowitsch     PetscCall(VecDot(auglag->Ciwork, auglag->Ciwork, &ineq_norm)); /* contribution to scalar Lagrangian */
656661095bbSAlp Dener     /* dL/dX += mu * pmax(0, Ci + Yi/mu)^T Ai */
6579566063dSJacob Faibussowitsch     PetscCall(VecScale(auglag->Ciwork, auglag->mu));
6589566063dSJacob Faibussowitsch     PetscCall(MatMultTransposeAdd(auglag->Ai, auglag->Ciwork, auglag->LgradX, auglag->LgradX));
659a5b23f4aSJose E. Roman     /* dL/dS = 0 because there are no slacks in PHR */
6609566063dSJacob Faibussowitsch     PetscCall(VecZeroEntries(auglag->LgradS));
661661095bbSAlp Dener   }
662661095bbSAlp Dener   /* combine gradient together */
6639566063dSJacob Faibussowitsch   PetscCall(TaoALMMCombinePrimal_Private(tao, auglag->LgradX, auglag->LgradS, auglag->G));
664661095bbSAlp 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)] */
6655f80ce2aSJacob Faibussowitsch   auglag->Lval = auglag->fval + 0.5*auglag->mu*(eq_norm + ineq_norm);
666661095bbSAlp Dener   PetscFunctionReturn(0);
667661095bbSAlp Dener }
668661095bbSAlp Dener 
669661095bbSAlp Dener /*
670661095bbSAlp Dener Lc = F + Ye^TCe + Yi^T(Ci - S) + 0.5*mu*[Ce^TCe + (Ci - S)^T(Ci - S)]
671661095bbSAlp Dener 
672661095bbSAlp Dener dLc/dX = dF/dX + Ye^TAe + Yi^TAi + 0.5*mu*[Ce^TAe + (Ci - S)^TAi]
673661095bbSAlp Dener 
674661095bbSAlp Dener dLc/dS = -[Yi + mu*(Ci - S)]
675661095bbSAlp Dener */
676661095bbSAlp Dener static PetscErrorCode TaoALMMComputeAugLagAndGradient_Private(Tao tao)
677661095bbSAlp Dener {
678661095bbSAlp Dener   TAO_ALMM       *auglag = (TAO_ALMM*)tao->data;
679661095bbSAlp Dener   PetscReal      yeTce=0.0, yiTcims=0.0, ceTce=0.0, cimsTcims=0.0;
680661095bbSAlp Dener 
681661095bbSAlp Dener   PetscFunctionBegin;
6829566063dSJacob Faibussowitsch   PetscCall(TaoALMMEvaluateIterate_Private(tao, auglag->P));
683661095bbSAlp Dener   if (tao->eq_constrained) {
684661095bbSAlp Dener     /* compute scalar contributions */
6859566063dSJacob Faibussowitsch     PetscCall(VecDot(auglag->Ye, auglag->Ce, &yeTce));
6869566063dSJacob Faibussowitsch     PetscCall(VecDot(auglag->Ce, auglag->Ce, &ceTce));
687661095bbSAlp Dener     /* dL/dX += ye^T Ae */
6889566063dSJacob Faibussowitsch     PetscCall(MatMultTransposeAdd(auglag->Ae, auglag->Ye, auglag->LgradX, auglag->LgradX));
689661095bbSAlp Dener     /* dL/dX += mu * ce^T Ae */
6909566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(auglag->Ae, auglag->Ce, auglag->Xwork));
6919566063dSJacob Faibussowitsch     PetscCall(VecAXPY(auglag->LgradX, auglag->mu, auglag->Xwork));
692661095bbSAlp Dener   }
693661095bbSAlp Dener   if (tao->ineq_constrained) {
694661095bbSAlp Dener     /* compute scalar contributions */
6959566063dSJacob Faibussowitsch     PetscCall(VecDot(auglag->Yi, auglag->Ci, &yiTcims));
6969566063dSJacob Faibussowitsch     PetscCall(VecDot(auglag->Ci, auglag->Ci, &cimsTcims));
697661095bbSAlp Dener     /* dL/dX += yi^T Ai */
6989566063dSJacob Faibussowitsch     PetscCall(MatMultTransposeAdd(auglag->Ai, auglag->Yi, auglag->LgradX, auglag->LgradX));
699661095bbSAlp Dener     /* dL/dX += mu * (ci - s)^T Ai */
7009566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(auglag->Ai, auglag->Ci, auglag->Xwork));
7019566063dSJacob Faibussowitsch     PetscCall(VecAXPY(auglag->LgradX, auglag->mu, auglag->Xwork));
702661095bbSAlp Dener     /* dL/dS = -[yi + mu*(ci - s)] */
7039566063dSJacob Faibussowitsch     PetscCall(VecWAXPY(auglag->LgradS, auglag->mu, auglag->Ci, auglag->Yi));
7049566063dSJacob Faibussowitsch     PetscCall(VecScale(auglag->LgradS, -1.0));
705661095bbSAlp Dener   }
706661095bbSAlp Dener   /* combine gradient together */
7079566063dSJacob Faibussowitsch   PetscCall(TaoALMMCombinePrimal_Private(tao, auglag->LgradX, auglag->LgradS, auglag->G));
708661095bbSAlp Dener   /* compute L = f + ye^T ce + yi^T (ci - s) + 0.5*mu*||ce||^2 + 0.5*mu*||ci - s||^2 */
709661095bbSAlp Dener   auglag->Lval = auglag->fval + yeTce + yiTcims + 0.5*auglag->mu*(ceTce + cimsTcims);
710661095bbSAlp Dener   PetscFunctionReturn(0);
711661095bbSAlp Dener }
712661095bbSAlp Dener 
7137721a36fSStefano Zampini PetscErrorCode TaoALMMSubsolverObjective_Private(Tao tao, Vec P, PetscReal *Lval, void *ctx)
7147721a36fSStefano Zampini {
7157721a36fSStefano Zampini   TAO_ALMM       *auglag = (TAO_ALMM*)ctx;
7167721a36fSStefano Zampini 
7177721a36fSStefano Zampini   PetscFunctionBegin;
7189566063dSJacob Faibussowitsch   PetscCall(VecCopy(P, auglag->P));
7199566063dSJacob Faibussowitsch   PetscCall((*auglag->sub_obj)(auglag->parent));
7207721a36fSStefano Zampini   *Lval = auglag->Lval;
7217721a36fSStefano Zampini   PetscFunctionReturn(0);
7227721a36fSStefano Zampini }
7237721a36fSStefano Zampini 
724661095bbSAlp Dener PetscErrorCode TaoALMMSubsolverObjectiveAndGradient_Private(Tao tao, Vec P, PetscReal *Lval, Vec G, void *ctx)
725661095bbSAlp Dener {
726661095bbSAlp Dener   TAO_ALMM       *auglag = (TAO_ALMM*)ctx;
727661095bbSAlp Dener 
728661095bbSAlp Dener   PetscFunctionBegin;
7299566063dSJacob Faibussowitsch   PetscCall(VecCopy(P, auglag->P));
7309566063dSJacob Faibussowitsch   PetscCall((*auglag->sub_obj)(auglag->parent));
7319566063dSJacob Faibussowitsch   PetscCall(VecCopy(auglag->G, G));
732661095bbSAlp Dener   *Lval = auglag->Lval;
733661095bbSAlp Dener   PetscFunctionReturn(0);
734661095bbSAlp Dener }
735