xref: /petsc/src/tao/constrained/impls/almm/almm.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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) {
23*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecZeroEntries(auglag->Ps));
24*5f80ce2aSJacob Faibussowitsch       CHKERRQ(TaoALMMCombinePrimal_Private(tao, auglag->Px, auglag->Ps, auglag->P));
25*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecZeroEntries(auglag->Yi));
26661095bbSAlp Dener     }
27661095bbSAlp Dener     if (tao->eq_constrained) {
28*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecZeroEntries(auglag->Ye));
29661095bbSAlp Dener     }
30661095bbSAlp Dener   }
31661095bbSAlp Dener 
32661095bbSAlp Dener   /* compute initial nonlinear Lagrangian and its derivatives */
33*5f80ce2aSJacob Faibussowitsch   CHKERRQ((*auglag->sub_obj)(tao));
34*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoALMMComputeOptimalityNorms_Private(tao));
35661095bbSAlp Dener   /* print initial step and check convergence */
36*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscInfo(tao,"Solving with %s formulation\n",TaoALMMTypes[auglag->type]));
37*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoLogConvergenceHistory(tao, auglag->Lval, auglag->gnorm, auglag->cnorm, tao->ksp_its));
38*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoMonitor(tao, tao->niter, auglag->fval, auglag->gnorm, auglag->cnorm, 0.0));
39*5f80ce2aSJacob Faibussowitsch   CHKERRQ((*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) {
48*5f80ce2aSJacob Faibussowitsch         CHKERRQ(VecDot(auglag->Ce, auglag->Ce, &auglag->cenorm));
49661095bbSAlp Dener       }
50661095bbSAlp Dener       auglag->cinorm = 0.0;
51661095bbSAlp Dener       if (tao->ineq_constrained) {
52*5f80ce2aSJacob Faibussowitsch         CHKERRQ(VecCopy(auglag->Ci, auglag->Ciwork));
53*5f80ce2aSJacob Faibussowitsch         CHKERRQ(VecScale(auglag->Ciwork, -1.0));
54*5f80ce2aSJacob Faibussowitsch         CHKERRQ(VecPointwiseMax(auglag->Ciwork, auglag->Cizero, auglag->Ciwork));
55*5f80ce2aSJacob Faibussowitsch         CHKERRQ(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;
64*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscInfo(tao,"Initial penalty: %.2f\n",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 */
70*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscInfo(tao,"Subsolver tolerance: ||G|| <= %e\n",auglag->gtol));
71*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TaoSetTolerances(auglag->subsolver, auglag->gtol, 0.0, 0.0));
72661095bbSAlp Dener     /* solve the bound-constrained or unconstrained subproblem */
73*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TaoSolve(auglag->subsolver));
74*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TaoGetConvergedReason(auglag->subsolver, &reason));
75661095bbSAlp Dener     tao->ksp_its += auglag->subsolver->ksp_its;
76661095bbSAlp Dener     if (reason != TAO_CONVERGED_GATOL) {
77*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscInfo(tao,"Subsolver failed to converge, reason: %s\n",TaoConvergedReasons[reason]));
78661095bbSAlp Dener     }
79661095bbSAlp Dener     /* evaluate solution and test convergence */
80*5f80ce2aSJacob Faibussowitsch     CHKERRQ((*auglag->sub_obj)(tao));
81*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TaoALMMComputeOptimalityNorms_Private(tao));
82661095bbSAlp Dener     /* decide whether to update multipliers or not */
83661095bbSAlp Dener     updated = 0.0;
84661095bbSAlp Dener     if (auglag->cnorm <= auglag->ytol) {
85*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscInfo(tao,"Multipliers updated: ||C|| <= %e\n",auglag->ytol));
86661095bbSAlp Dener       /* constraints are good, update multipliers and convergence tolerances */
87661095bbSAlp Dener       if (tao->eq_constrained) {
88*5f80ce2aSJacob Faibussowitsch         CHKERRQ(VecAXPY(auglag->Ye, auglag->mu, auglag->Ce));
89*5f80ce2aSJacob Faibussowitsch         CHKERRQ(VecSet(auglag->Cework, auglag->ye_max));
90*5f80ce2aSJacob Faibussowitsch         CHKERRQ(VecPointwiseMin(auglag->Ye, auglag->Cework, auglag->Ye));
91*5f80ce2aSJacob Faibussowitsch         CHKERRQ(VecSet(auglag->Cework, auglag->ye_min));
92*5f80ce2aSJacob Faibussowitsch         CHKERRQ(VecPointwiseMax(auglag->Ye, auglag->Cework, auglag->Ye));
93661095bbSAlp Dener       }
94661095bbSAlp Dener       if (tao->ineq_constrained) {
95*5f80ce2aSJacob Faibussowitsch         CHKERRQ(VecAXPY(auglag->Yi, auglag->mu, auglag->Ci));
96*5f80ce2aSJacob Faibussowitsch         CHKERRQ(VecSet(auglag->Ciwork, auglag->yi_max));
97*5f80ce2aSJacob Faibussowitsch         CHKERRQ(VecPointwiseMin(auglag->Yi, auglag->Ciwork, auglag->Yi));
98*5f80ce2aSJacob Faibussowitsch         CHKERRQ(VecSet(auglag->Ciwork, auglag->yi_min));
99*5f80ce2aSJacob Faibussowitsch         CHKERRQ(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       }
115*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscInfo(tao,"Penalty increased: mu = %.2f\n",auglag->mu));
116661095bbSAlp Dener     }
117*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TaoLogConvergenceHistory(tao, auglag->fval, auglag->gnorm, auglag->cnorm, tao->ksp_its));
118*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TaoMonitor(tao, tao->niter, auglag->fval, auglag->gnorm, auglag->cnorm, updated));
119*5f80ce2aSJacob Faibussowitsch     CHKERRQ((*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;
131*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPushTab(viewer));
132*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoView(auglag->subsolver,viewer));
133*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerASCIIPopTab(viewer));
134*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
135661095bbSAlp Dener   if (isascii) {
136*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPushTab(viewer));
137*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(viewer, "ALMM Formulation Type: %s\n", TaoALMMTypes[auglag->type]));
138*5f80ce2aSJacob Faibussowitsch     CHKERRQ(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.");
153*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoComputeVariableBounds(tao));
154661095bbSAlp Dener   /* alias base vectors and create extras */
155*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetType(tao->solution, &vec_type));
156661095bbSAlp Dener   auglag->Px = tao->solution;
157661095bbSAlp Dener   if (!tao->gradient) { /* base gradient */
158*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDuplicate(tao->solution, &tao->gradient));
159661095bbSAlp Dener   }
160661095bbSAlp Dener   auglag->LgradX = tao->gradient;
161661095bbSAlp Dener   if (!auglag->Xwork) { /* opt var work vector */
162*5f80ce2aSJacob Faibussowitsch     CHKERRQ(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 */
168*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecDuplicate(auglag->Ce, &auglag->Ye));
169661095bbSAlp Dener     }
170661095bbSAlp Dener     if (!auglag->Cework) {
171*5f80ce2aSJacob Faibussowitsch       CHKERRQ(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 */
178*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecDuplicate(auglag->Ci, &auglag->Yi));
179661095bbSAlp Dener     }
180661095bbSAlp Dener     if (!auglag->Ciwork) {
181*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecDuplicate(auglag->Ci, &auglag->Ciwork));
182661095bbSAlp Dener     }
183661095bbSAlp Dener     if (!auglag->Cizero) {
184*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecDuplicate(auglag->Ci, &auglag->Cizero));
185*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecZeroEntries(auglag->Cizero));
186661095bbSAlp Dener     }
187661095bbSAlp Dener     if (!auglag->Ps) { /* slack vars */
188*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecDuplicate(auglag->Ci, &auglag->Ps));
189661095bbSAlp Dener     }
190661095bbSAlp Dener     if (!auglag->LgradS) { /* slack component of Lagrangian gradient */
191*5f80ce2aSJacob Faibussowitsch       CHKERRQ(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) {
195*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscMalloc1(2, &auglag->Parr));
196661095bbSAlp Dener       auglag->Parr[0] = auglag->Px; auglag->Parr[1] = auglag->Ps;
197*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecConcatenate(2, auglag->Parr, &auglag->P, &auglag->Pis));
198*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscMalloc1(2, &auglag->Pscatter));
199*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecScatterCreate(auglag->P, auglag->Pis[0], auglag->Px, NULL, &auglag->Pscatter[0]));
200*5f80ce2aSJacob Faibussowitsch       CHKERRQ(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) {
205*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscMalloc1(2, &auglag->Yarr));
206661095bbSAlp Dener         auglag->Yarr[0] = auglag->Ye; auglag->Yarr[1] = auglag->Yi;
207*5f80ce2aSJacob Faibussowitsch         CHKERRQ(VecConcatenate(2, auglag->Yarr, &auglag->Y, &auglag->Yis));
208*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscMalloc1(2, &auglag->Yscatter));
209*5f80ce2aSJacob Faibussowitsch         CHKERRQ(VecScatterCreate(auglag->Y, auglag->Yis[0], auglag->Ye, NULL, &auglag->Yscatter[0]));
210*5f80ce2aSJacob Faibussowitsch         CHKERRQ(VecScatterCreate(auglag->Y, auglag->Yis[1], auglag->Yi, NULL, &auglag->Yscatter[1]));
211661095bbSAlp Dener       }
212661095bbSAlp Dener       if (!auglag->C) {
213*5f80ce2aSJacob Faibussowitsch         CHKERRQ(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) {
244*5f80ce2aSJacob Faibussowitsch       CHKERRQ(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 */
260*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetSolution(auglag->subsolver, auglag->P));
261*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetObjective(auglag->subsolver, TaoALMMSubsolverObjective_Private, (void*)auglag));
262*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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 */
265*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectTypeCompare((PetscObject)auglag->subsolver, TAOCG, &is_cg));
266*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectTypeCompare((PetscObject)auglag->subsolver, TAOLMVM, &is_lmvm));
267661095bbSAlp Dener     if (is_cg) {
268*5f80ce2aSJacob Faibussowitsch       CHKERRQ(TaoSetType(auglag->subsolver, TAOBNCG));
269*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscInfo(tao,"TAOCG detected for bound-constrained problem, switching to TAOBNCG instead."));
270661095bbSAlp Dener     }
271661095bbSAlp Dener     if (is_lmvm) {
272*5f80ce2aSJacob Faibussowitsch       CHKERRQ(TaoSetType(auglag->subsolver, TAOBQNLS));
273*5f80ce2aSJacob Faibussowitsch       CHKERRQ(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) {
277*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecDuplicate(auglag->P, &auglag->PL));
278661095bbSAlp Dener     }
279661095bbSAlp Dener     if (!auglag->PU) {
280*5f80ce2aSJacob Faibussowitsch       CHKERRQ(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 */
284*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecDuplicate(auglag->Ci, &SL));
285*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecSet(SL, 0.0));
286*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecDuplicate(auglag->Ci, &SU));
287*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecSet(SU, PETSC_INFINITY));
288661095bbSAlp Dener       /* combine opt var bounds with slack bounds */
289*5f80ce2aSJacob Faibussowitsch       CHKERRQ(TaoALMMCombinePrimal_Private(tao, tao->XL, SL, auglag->PL));
290*5f80ce2aSJacob Faibussowitsch       CHKERRQ(TaoALMMCombinePrimal_Private(tao, tao->XU, SU, auglag->PU));
291661095bbSAlp Dener       /* destroy work vectors */
292*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecDestroy(&SL));
293*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecDestroy(&SU));
294661095bbSAlp Dener     } else {
295661095bbSAlp Dener       /* no inequality constraints, just copy bounds into the subsolver */
296*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecCopy(tao->XL, auglag->PL));
297*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecCopy(tao->XU, auglag->PU));
298661095bbSAlp Dener     }
299*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TaoSetVariableBounds(auglag->subsolver, auglag->PL, auglag->PU));
300661095bbSAlp Dener   }
301*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
312*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoDestroy(&auglag->subsolver));
313661095bbSAlp Dener   if (tao->setupcalled) {
314*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&auglag->Xwork));              /* opt work */
315661095bbSAlp Dener     if (tao->eq_constrained) {
316*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecDestroy(&auglag->Ye));               /* equality multipliers */
317*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecDestroy(&auglag->Cework));           /* equality work vector */
318661095bbSAlp Dener     }
319661095bbSAlp Dener     if (tao->ineq_constrained) {
320*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecDestroy(&auglag->Ps));               /* slack vars */
321661095bbSAlp Dener       auglag->Parr[0] = NULL;                                     /* clear pointer to tao->solution, will be destroyed by TaoDestroy() shell */
322*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFree(auglag->Parr));               /* array of primal vectors */
323*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecDestroy(&auglag->LgradS));           /* slack grad */
324*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecDestroy(&auglag->Cizero));           /* zero vector for pointwise max */
325*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecDestroy(&auglag->Yi));               /* inequality multipliers */
326*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecDestroy(&auglag->Ciwork));           /* inequality work vector */
327*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecDestroy(&auglag->P));                /* combo primal */
328*5f80ce2aSJacob Faibussowitsch       CHKERRQ(ISDestroy(&auglag->Pis[0]));            /* index set for X inside P */
329*5f80ce2aSJacob Faibussowitsch       CHKERRQ(ISDestroy(&auglag->Pis[1]));            /* index set for S inside P */
330*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFree(auglag->Pis));                /* array of P index sets */
331*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecScatterDestroy(&auglag->Pscatter[0]));
332*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecScatterDestroy(&auglag->Pscatter[1]));
333*5f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFree(auglag->Pscatter));
334661095bbSAlp Dener       if (tao->eq_constrained) {
335*5f80ce2aSJacob Faibussowitsch         CHKERRQ(VecDestroy(&auglag->Y));              /* combo multipliers */
336*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscFree(auglag->Yarr));             /* array of dual vectors */
337*5f80ce2aSJacob Faibussowitsch         CHKERRQ(VecDestroy(&auglag->C));              /* combo constraints */
338*5f80ce2aSJacob Faibussowitsch         CHKERRQ(ISDestroy(&auglag->Yis[0]));          /* index set for Ye inside Y */
339*5f80ce2aSJacob Faibussowitsch         CHKERRQ(ISDestroy(&auglag->Yis[1]));          /* index set for Yi inside Y */
340*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscFree(auglag->Yis));
341*5f80ce2aSJacob Faibussowitsch         CHKERRQ(VecScatterDestroy(&auglag->Yscatter[0]));
342*5f80ce2aSJacob Faibussowitsch         CHKERRQ(VecScatterDestroy(&auglag->Yscatter[1]));
343*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscFree(auglag->Yscatter));
344661095bbSAlp Dener       }
345661095bbSAlp Dener     }
346661095bbSAlp Dener     if (tao->bounded) {
347*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecDestroy(&auglag->PL));                /* lower bounds for subsolver */
348*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecDestroy(&auglag->PU));                /* upper bounds for subsolver */
349661095bbSAlp Dener     }
350661095bbSAlp Dener   }
351*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(tao->data));
352661095bbSAlp Dener   PetscFunctionReturn(0);
353661095bbSAlp Dener }
354661095bbSAlp Dener 
355661095bbSAlp Dener static PetscErrorCode TaoSetFromOptions_ALMM(PetscOptionItems *PetscOptionsObject,Tao tao)
356661095bbSAlp Dener {
357661095bbSAlp Dener   TAO_ALMM       *auglag = (TAO_ALMM*)tao->data;
358661095bbSAlp Dener   PetscInt       i;
359661095bbSAlp Dener 
360661095bbSAlp Dener   PetscFunctionBegin;
361*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHead(PetscOptionsObject,"Augmented Lagrangian multipler method solves problems with general constraints by converting them into a sequence of unconstrained problems."));
362*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsReal("-tao_almm_mu_init","initial penalty parameter","",auglag->mu0,&auglag->mu0,NULL));
363*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsReal("-tao_almm_mu_factor","increase factor for the penalty parameter","",auglag->mu_fac,&auglag->mu_fac,NULL));
364*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsReal("-tao_almm_mu_power_good","exponential for penalty parameter when multiplier update is accepted","",auglag->mu_pow_good,&auglag->mu_pow_good,NULL));
365*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsReal("-tao_almm_mu_power_bad","exponential for penalty parameter when multiplier update is rejected","",auglag->mu_pow_bad,&auglag->mu_pow_bad,NULL));
366*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsReal("-tao_almm_mu_max","maximum safeguard for penalty parameter updates","",auglag->mu_max,&auglag->mu_max,NULL));
367*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsReal("-tao_almm_ye_min","minimum safeguard for equality multiplier updates","",auglag->ye_min,&auglag->ye_min,NULL));
368*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsReal("-tao_almm_ye_max","maximum safeguard for equality multipliers updates","",auglag->ye_max,&auglag->ye_max,NULL));
369*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsReal("-tao_almm_yi_min","minimum safeguard for inequality multipliers updates","",auglag->yi_min,&auglag->yi_min,NULL));
370*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsReal("-tao_almm_yi_max","maximum safeguard for inequality multipliers updates","",auglag->yi_max,&auglag->yi_max,NULL));
371*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsEnum("-tao_almm_type","augmented Lagrangian formulation type for the subproblem","TaoALMMType",TaoALMMTypes,(PetscEnum)auglag->type,(PetscEnum*)&auglag->type,NULL));
372*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsTail());
373*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetOptionsPrefix(auglag->subsolver,((PetscObject)tao)->prefix));
374*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoAppendOptionsPrefix(auglag->subsolver,"tao_almm_subsolver_"));
375*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetFromOptions(auglag->subsolver));
376661095bbSAlp Dener   for (i=0; i<tao->numbermonitors; i++) {
377*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectReference((PetscObject)tao->monitorcontext[i]));
378*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TaoSetMonitor(auglag->subsolver, tao->monitor[i], tao->monitorcontext[i], tao->monitordestroy[i]));
379661095bbSAlp Dener     if (tao->monitor[i] == TaoMonitorDefault || tao->monitor[i] == TaoDefaultCMonitor || tao->monitor[i] == TaoDefaultGMonitor || tao->monitor[i] == TaoDefaultSMonitor) {
380661095bbSAlp Dener       auglag->info = PETSC_TRUE;
381661095bbSAlp Dener     }
382661095bbSAlp Dener   }
383661095bbSAlp Dener   PetscFunctionReturn(0);
384661095bbSAlp Dener }
385661095bbSAlp Dener 
386661095bbSAlp Dener /* -------------------------------------------------------- */
387661095bbSAlp Dener 
388661095bbSAlp Dener /*MC
389661095bbSAlp Dener   TaoALMM - Augmented Lagrangian multiplier method for solving nonlinear optimization problems with general constraints.
390661095bbSAlp Dener 
391661095bbSAlp Dener   Options Database Keys:
392661095bbSAlp Dener + -tao_almm_mu_init <real>       - initial penalty parameter (default: 10.)
393661095bbSAlp Dener . -tao_almm_mu_factor <real>     - increase factor for the penalty parameter (default: 100.)
394661095bbSAlp Dener . -tao_almm_mu_max <real>        - maximum safeguard for penalty parameter updates (default: 1.e20)
395661095bbSAlp Dener . -tao_almm_mu_power_good <real> - exponential for penalty parameter when multiplier update is accepted (default: 0.9)
396661095bbSAlp Dener . -tao_almm_mu_power_bad <real>  - exponential for penalty parameter when multiplier update is rejected (default: 0.1)
397661095bbSAlp Dener . -tao_almm_ye_min <real>        - minimum safeguard for equality multiplier updates (default: -1.e20)
398661095bbSAlp Dener . -tao_almm_ye_max <real>        - maximum safeguard for equality multiplier updates (default: 1.e20)
399661095bbSAlp Dener . -tao_almm_yi_min <real>        - minimum safeguard for inequality multiplier updates (default: -1.e20)
400661095bbSAlp Dener . -tao_almm_yi_max <real>        - maximum safeguard for inequality multiplier updates (default: 1.e20)
401661095bbSAlp Dener - -tao_almm_type <classic,phr>   - change formulation of the augmented Lagrangian merit function for the subproblem (default: classic)
402661095bbSAlp Dener 
403661095bbSAlp Dener   Level: beginner
404661095bbSAlp Dener 
405661095bbSAlp Dener   Notes:
406661095bbSAlp Dener   This method converts a constrained problem into a sequence of unconstrained problems via the augmented
407661095bbSAlp Dener   Lagrangian merit function. Bound constraints are pushed down to the subproblem without any modifications.
408661095bbSAlp Dener 
409661095bbSAlp Dener   Two formulations are offered for the subproblem: canonical Hestenes-Powell augmented Lagrangian with slack
410661095bbSAlp Dener   variables for inequality constraints, and a slack-less Powell-Hestenes-Rockafellar (PHR) formulation utilizing a
411661095bbSAlp Dener   pointwise max() penalty on inequality constraints. The canonical augmented Lagrangian formulation typically
412661095bbSAlp Dener   converges faster for most problems. However, PHR may be desirable for problems featuring a large number
413661095bbSAlp Dener   of inequality constraints because it avoids inflating the size of the subproblem with slack variables.
414661095bbSAlp Dener 
415661095bbSAlp Dener   The subproblem is solved using a nested first-order TAO solver. The user can retrieve a pointer to
416661095bbSAlp Dener   the subsolver via TaoALMMGetSubsolver() or pass command line arguments to it using the
417661095bbSAlp Dener   "-tao_almm_subsolver_" prefix. Currently, TaoALMM does not support second-order methods for the
418661095bbSAlp Dener   subproblem. It is also highly recommended that the subsolver chosen by the user utilize a trust-region
419661095bbSAlp Dener   strategy for globalization (default: TAOBQNKTR) especially if the outer problem features bound constraints.
420661095bbSAlp Dener 
421661095bbSAlp Dener .vb
422661095bbSAlp Dener   while unconverged
423661095bbSAlp Dener     solve argmin_x L(x) s.t. l <= x <= u
424661095bbSAlp Dener     if ||c|| <= y_tol
425661095bbSAlp Dener       if ||c|| <= c_tol && ||Lgrad|| <= g_tol:
426661095bbSAlp Dener         problem converged, return solution
427661095bbSAlp Dener       else
428661095bbSAlp Dener         constraints sufficiently improved
429661095bbSAlp Dener         update multipliers and tighten tolerances
430661095bbSAlp Dener       endif
431661095bbSAlp Dener     else
432661095bbSAlp Dener       constraints did not improve
433661095bbSAlp Dener       update penalty and loosen tolerances
434661095bbSAlp Dener     endif
435661095bbSAlp Dener   endwhile
436661095bbSAlp Dener .ve
437661095bbSAlp Dener 
438661095bbSAlp Dener .seealso: TaoALMMGetType(), TaoALMMSetType(), TaoALMMSetSubsolver(), TaoALMMGetSubsolver(),
439661095bbSAlp Dener           TaoALMMGetMultipliers(), TaoALMMSetMultipliers(), TaoALMMGetPrimalIS(), TaoALMMGetDualIS()
440661095bbSAlp Dener M*/
441661095bbSAlp Dener PETSC_EXTERN PetscErrorCode TaoCreate_ALMM(Tao tao)
442661095bbSAlp Dener {
443661095bbSAlp Dener   TAO_ALMM       *auglag;
444661095bbSAlp Dener 
445661095bbSAlp Dener   PetscFunctionBegin;
446*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscNewLog(tao, &auglag));
447661095bbSAlp Dener 
448661095bbSAlp Dener   tao->ops->destroy        = TaoDestroy_ALMM;
449661095bbSAlp Dener   tao->ops->setup          = TaoSetUp_ALMM;
450661095bbSAlp Dener   tao->ops->setfromoptions = TaoSetFromOptions_ALMM;
451661095bbSAlp Dener   tao->ops->view           = TaoView_ALMM;
452661095bbSAlp Dener   tao->ops->solve          = TaoSolve_ALMM;
453661095bbSAlp Dener 
454661095bbSAlp Dener   tao->gatol = 1.e-5;
455661095bbSAlp Dener   tao->grtol = 0.0;
456661095bbSAlp Dener   tao->gttol = 0.0;
457661095bbSAlp Dener   tao->catol = 1.e-5;
458661095bbSAlp Dener   tao->crtol = 0.0;
459661095bbSAlp Dener 
460661095bbSAlp Dener   tao->data           = (void*)auglag;
461661095bbSAlp Dener   auglag->parent      = tao;
462661095bbSAlp Dener   auglag->mu0         = 10.0;
463661095bbSAlp Dener   auglag->mu          = auglag->mu0;
464661095bbSAlp Dener   auglag->mu_fac      = 10.0;
465661095bbSAlp Dener   auglag->mu_max      = PETSC_INFINITY;
466661095bbSAlp Dener   auglag->mu_pow_good = 0.9;
467661095bbSAlp Dener   auglag->mu_pow_bad  = 0.1;
468661095bbSAlp Dener   auglag->ye_min      = PETSC_NINFINITY;
469661095bbSAlp Dener   auglag->ye_max      = PETSC_INFINITY;
470661095bbSAlp Dener   auglag->yi_min      = PETSC_NINFINITY;
471661095bbSAlp Dener   auglag->yi_max      = PETSC_INFINITY;
472661095bbSAlp Dener   auglag->ytol0       = 0.1/PetscPowReal(auglag->mu0, auglag->mu_pow_bad);
473661095bbSAlp Dener   auglag->ytol        = auglag->ytol0;
474661095bbSAlp Dener   auglag->gtol0       = 1.0/auglag->mu0;
475661095bbSAlp Dener   auglag->gtol        = auglag->gtol0;
476661095bbSAlp Dener 
477661095bbSAlp Dener   auglag->sub_obj     = TaoALMMComputeAugLagAndGradient_Private;
478661095bbSAlp Dener   auglag->type        = TAO_ALMM_CLASSIC;
479661095bbSAlp Dener   auglag->info        = PETSC_FALSE;
480661095bbSAlp Dener 
481*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoCreate(PetscObjectComm((PetscObject)tao),&auglag->subsolver));
482*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetType(auglag->subsolver, TAOBQNKTR));
483*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetTolerances(auglag->subsolver, auglag->gtol, 0.0, 0.0));
484*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetMaximumIterations(auglag->subsolver, 1000));
485*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetMaximumFunctionEvaluations(auglag->subsolver, 10000));
486*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetFunctionLowerBound(auglag->subsolver, PETSC_NINFINITY));
487*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectIncrementTabLevel((PetscObject)auglag->subsolver,(PetscObject)tao,1));
488661095bbSAlp Dener 
489*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetType_C", TaoALMMGetType_Private));
490*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMSetType_C", TaoALMMSetType_Private));
491*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetSubsolver_C", TaoALMMGetSubsolver_Private));
492*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMSetSubsolver_C", TaoALMMSetSubsolver_Private));
493*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetMultipliers_C", TaoALMMGetMultipliers_Private));
494*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMSetMultipliers_C", TaoALMMSetMultipliers_Private));
495*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetPrimalIS_C", TaoALMMGetPrimalIS_Private));
496*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetDualIS_C", TaoALMMGetDualIS_Private));
497661095bbSAlp Dener   PetscFunctionReturn(0);
498661095bbSAlp Dener }
499661095bbSAlp Dener 
500661095bbSAlp Dener static PetscErrorCode TaoALMMCombinePrimal_Private(Tao tao, Vec X, Vec S, Vec P)
501661095bbSAlp Dener {
502661095bbSAlp Dener   TAO_ALMM     *auglag = (TAO_ALMM*)tao->data;
503661095bbSAlp Dener 
504661095bbSAlp Dener   PetscFunctionBegin;
505661095bbSAlp Dener   if (tao->ineq_constrained) {
506*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterBegin(auglag->Pscatter[0], X, P, INSERT_VALUES, SCATTER_REVERSE));
507*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterEnd(auglag->Pscatter[0], X, P, INSERT_VALUES, SCATTER_REVERSE));
508*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterBegin(auglag->Pscatter[1], S, P, INSERT_VALUES, SCATTER_REVERSE));
509*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterEnd(auglag->Pscatter[1], S, P, INSERT_VALUES, SCATTER_REVERSE));
510661095bbSAlp Dener   } else {
511*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCopy(X, P));
512661095bbSAlp Dener   }
513661095bbSAlp Dener   PetscFunctionReturn(0);
514661095bbSAlp Dener }
515661095bbSAlp Dener 
516661095bbSAlp Dener static PetscErrorCode TaoALMMCombineDual_Private(Tao tao, Vec EQ, Vec IN, Vec Y)
517661095bbSAlp Dener {
518661095bbSAlp Dener   TAO_ALMM     *auglag = (TAO_ALMM*)tao->data;
519661095bbSAlp Dener 
520661095bbSAlp Dener   PetscFunctionBegin;
521661095bbSAlp Dener   if (tao->eq_constrained) {
522661095bbSAlp Dener     if (tao->ineq_constrained) {
523*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecScatterBegin(auglag->Yscatter[0], EQ, Y, INSERT_VALUES, SCATTER_REVERSE));
524*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecScatterEnd(auglag->Yscatter[0], EQ, Y, INSERT_VALUES, SCATTER_REVERSE));
525*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecScatterBegin(auglag->Yscatter[1], IN, Y, INSERT_VALUES, SCATTER_REVERSE));
526*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecScatterEnd(auglag->Yscatter[1], IN, Y, INSERT_VALUES, SCATTER_REVERSE));
527661095bbSAlp Dener     } else {
528*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecCopy(EQ, Y));
529661095bbSAlp Dener     }
530661095bbSAlp Dener   } else {
531*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCopy(IN, Y));
532661095bbSAlp Dener   }
533661095bbSAlp Dener   PetscFunctionReturn(0);
534661095bbSAlp Dener }
535661095bbSAlp Dener 
536661095bbSAlp Dener static PetscErrorCode TaoALMMSplitPrimal_Private(Tao tao, Vec P, Vec X, Vec S)
537661095bbSAlp Dener {
538661095bbSAlp Dener   TAO_ALMM     *auglag = (TAO_ALMM*)tao->data;
539661095bbSAlp Dener 
540661095bbSAlp Dener   PetscFunctionBegin;
541661095bbSAlp Dener   if (tao->ineq_constrained) {
542*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterBegin(auglag->Pscatter[0], P, X, INSERT_VALUES, SCATTER_FORWARD));
543*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterEnd(auglag->Pscatter[0], P, X, INSERT_VALUES, SCATTER_FORWARD));
544*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterBegin(auglag->Pscatter[1], P, S, INSERT_VALUES, SCATTER_FORWARD));
545*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScatterEnd(auglag->Pscatter[1], P, S, INSERT_VALUES, SCATTER_FORWARD));
546661095bbSAlp Dener   } else {
547*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCopy(P, X));
548661095bbSAlp Dener   }
549661095bbSAlp Dener   PetscFunctionReturn(0);
550661095bbSAlp Dener }
551661095bbSAlp Dener 
552661095bbSAlp Dener /* this assumes that the latest constraints are stored in Ce and Ci, and also combined in C */
553661095bbSAlp Dener static PetscErrorCode TaoALMMComputeOptimalityNorms_Private(Tao tao)
554661095bbSAlp Dener {
555661095bbSAlp Dener   TAO_ALMM     *auglag = (TAO_ALMM*)tao->data;
556661095bbSAlp Dener 
557661095bbSAlp Dener   PetscFunctionBegin;
558661095bbSAlp Dener   /* if bounded, project the gradient */
559661095bbSAlp Dener   if (tao->bounded) {
560*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecBoundGradientProjection(auglag->LgradX, auglag->Px, tao->XL, tao->XU, auglag->LgradX));
561661095bbSAlp Dener   }
562661095bbSAlp Dener   if (auglag->type == TAO_ALMM_PHR) {
563*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecNorm(auglag->LgradX, NORM_INFINITY, &auglag->gnorm));
564661095bbSAlp Dener     auglag->cenorm = 0.0;
565661095bbSAlp Dener     if (tao->eq_constrained) {
566*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecNorm(auglag->Ce, NORM_INFINITY, &auglag->cenorm));
567661095bbSAlp Dener     }
568661095bbSAlp Dener     auglag->cinorm = 0.0;
569661095bbSAlp Dener     if (tao->ineq_constrained) {
570*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecCopy(auglag->Yi, auglag->Ciwork));
571*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecScale(auglag->Ciwork, -1.0/auglag->mu));
572*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecPointwiseMax(auglag->Ciwork, auglag->Ci, auglag->Ciwork));
573*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecNorm(auglag->Ciwork, NORM_INFINITY, &auglag->cinorm));
574661095bbSAlp Dener     }
575661095bbSAlp Dener     auglag->cnorm_old = auglag->cnorm;
576661095bbSAlp Dener     auglag->cnorm = PetscMax(auglag->cenorm, auglag->cinorm);
577661095bbSAlp Dener     auglag->ytol = auglag->ytol0 * auglag->cnorm_old;
578661095bbSAlp Dener   } else {
579*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecNorm(auglag->LgradX, NORM_2, &auglag->gnorm));
580*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecNorm(auglag->C, NORM_2, &auglag->cnorm));
581661095bbSAlp Dener   }
582661095bbSAlp Dener   PetscFunctionReturn(0);
583661095bbSAlp Dener }
584661095bbSAlp Dener 
585661095bbSAlp Dener static PetscErrorCode TaoALMMEvaluateIterate_Private(Tao tao, Vec P)
586661095bbSAlp Dener {
587661095bbSAlp Dener   TAO_ALMM       *auglag = (TAO_ALMM*)tao->data;
588661095bbSAlp Dener 
589661095bbSAlp Dener   PetscFunctionBegin;
590661095bbSAlp Dener   /* split solution into primal and slack components */
591*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoALMMSplitPrimal_Private(tao, auglag->P, auglag->Px, auglag->Ps));
592661095bbSAlp Dener 
593661095bbSAlp Dener   /* compute f, df/dx and the constraints */
594*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoComputeObjectiveAndGradient(tao, auglag->Px, &auglag->fval, auglag->LgradX));
595661095bbSAlp Dener   if (tao->eq_constrained) {
596*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TaoComputeEqualityConstraints(tao, auglag->Px, auglag->Ce));
597*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TaoComputeJacobianEquality(tao, auglag->Px, auglag->Ae, auglag->Ae));
598661095bbSAlp Dener   }
599661095bbSAlp Dener   if (tao->ineq_constrained) {
600*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TaoComputeInequalityConstraints(tao, auglag->Px, auglag->Ci));
601*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TaoComputeJacobianInequality(tao, auglag->Px, auglag->Ai, auglag->Ai));
602661095bbSAlp Dener     switch (auglag->type) {
6037721a36fSStefano Zampini       case TAO_ALMM_CLASSIC:
604661095bbSAlp Dener         /* classic formulation converts inequality to equality constraints via slack variables */
605*5f80ce2aSJacob Faibussowitsch         CHKERRQ(VecAXPY(auglag->Ci, -1.0, auglag->Ps));
606661095bbSAlp Dener         break;
6077721a36fSStefano Zampini       case TAO_ALMM_PHR:
608661095bbSAlp Dener         /* PHR is based on Ci <= 0 while TAO defines Ci >= 0 so we hit it with a negative sign */
609*5f80ce2aSJacob Faibussowitsch         CHKERRQ(VecScale(auglag->Ci, -1.0));
610*5f80ce2aSJacob Faibussowitsch         CHKERRQ(MatScale(auglag->Ai, -1.0));
611661095bbSAlp Dener         break;
612661095bbSAlp Dener       default:
613661095bbSAlp Dener         break;
614661095bbSAlp Dener     }
615661095bbSAlp Dener   }
616661095bbSAlp Dener   /* combine constraints into one vector */
617*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoALMMCombineDual_Private(tao, auglag->Ce, auglag->Ci, auglag->C));
618661095bbSAlp Dener   PetscFunctionReturn(0);
619661095bbSAlp Dener }
620661095bbSAlp Dener 
621661095bbSAlp Dener /*
622661095bbSAlp Dener Lphr = f + 0.5*mu*[ (Ce + Ye/mu)^T (Ce + Ye/mu) + pmin(0, Ci + Yi/mu)^T pmin(0, Ci + Yi/mu)]
623661095bbSAlp Dener 
624661095bbSAlp Dener dLphr/dX = dF/dX + mu*[ (Ce + Ye/mu)^T Ae + pmin(0, Ci + Yi/mu)^T Ai]
625661095bbSAlp Dener 
626661095bbSAlp Dener dLphr/dS = 0
627661095bbSAlp Dener */
628661095bbSAlp Dener static PetscErrorCode TaoALMMComputePHRLagAndGradient_Private(Tao tao)
629661095bbSAlp Dener {
630661095bbSAlp Dener   TAO_ALMM       *auglag = (TAO_ALMM*)tao->data;
631661095bbSAlp Dener   PetscReal      eq_norm=0.0, ineq_norm=0.0;
632661095bbSAlp Dener 
633661095bbSAlp Dener   PetscFunctionBegin;
634*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoALMMEvaluateIterate_Private(tao, auglag->P));
635661095bbSAlp Dener   if (tao->eq_constrained) {
636661095bbSAlp Dener     /* Ce_work = mu*(Ce + Ye/mu) */
637*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecWAXPY(auglag->Cework, 1.0/auglag->mu, auglag->Ye, auglag->Ce));
638*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDot(auglag->Cework, auglag->Cework, &eq_norm)); /* contribution to scalar Lagrangian */
639*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScale(auglag->Cework, auglag->mu));
640661095bbSAlp Dener     /* dL/dX += mu*(Ce + Ye/mu)^T Ae */
641*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMultTransposeAdd(auglag->Ae, auglag->Cework, auglag->LgradX, auglag->LgradX));
642661095bbSAlp Dener   }
643661095bbSAlp Dener   if (tao->ineq_constrained) {
644661095bbSAlp Dener     /* Ci_work = mu * pmax(0, Ci + Yi/mu) where pmax() is pointwise max() */
645*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecWAXPY(auglag->Ciwork, 1.0/auglag->mu, auglag->Yi, auglag->Ci));
646*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecPointwiseMax(auglag->Ciwork, auglag->Cizero, auglag->Ciwork));
647*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDot(auglag->Ciwork, auglag->Ciwork, &ineq_norm)); /* contribution to scalar Lagrangian */
648661095bbSAlp Dener     /* dL/dX += mu * pmax(0, Ci + Yi/mu)^T Ai */
649*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScale(auglag->Ciwork, auglag->mu));
650*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMultTransposeAdd(auglag->Ai, auglag->Ciwork, auglag->LgradX, auglag->LgradX));
651a5b23f4aSJose E. Roman     /* dL/dS = 0 because there are no slacks in PHR */
652*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecZeroEntries(auglag->LgradS));
653661095bbSAlp Dener   }
654661095bbSAlp Dener   /* combine gradient together */
655*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoALMMCombinePrimal_Private(tao, auglag->LgradX, auglag->LgradS, auglag->G));
656661095bbSAlp 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)] */
657*5f80ce2aSJacob Faibussowitsch   auglag->Lval = auglag->fval + 0.5*auglag->mu*(eq_norm + ineq_norm);
658661095bbSAlp Dener   PetscFunctionReturn(0);
659661095bbSAlp Dener }
660661095bbSAlp Dener 
661661095bbSAlp Dener /*
662661095bbSAlp Dener Lc = F + Ye^TCe + Yi^T(Ci - S) + 0.5*mu*[Ce^TCe + (Ci - S)^T(Ci - S)]
663661095bbSAlp Dener 
664661095bbSAlp Dener dLc/dX = dF/dX + Ye^TAe + Yi^TAi + 0.5*mu*[Ce^TAe + (Ci - S)^TAi]
665661095bbSAlp Dener 
666661095bbSAlp Dener dLc/dS = -[Yi + mu*(Ci - S)]
667661095bbSAlp Dener */
668661095bbSAlp Dener static PetscErrorCode TaoALMMComputeAugLagAndGradient_Private(Tao tao)
669661095bbSAlp Dener {
670661095bbSAlp Dener   TAO_ALMM       *auglag = (TAO_ALMM*)tao->data;
671661095bbSAlp Dener   PetscReal      yeTce=0.0, yiTcims=0.0, ceTce=0.0, cimsTcims=0.0;
672661095bbSAlp Dener 
673661095bbSAlp Dener   PetscFunctionBegin;
674*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoALMMEvaluateIterate_Private(tao, auglag->P));
675661095bbSAlp Dener   if (tao->eq_constrained) {
676661095bbSAlp Dener     /* compute scalar contributions */
677*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDot(auglag->Ye, auglag->Ce, &yeTce));
678*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDot(auglag->Ce, auglag->Ce, &ceTce));
679661095bbSAlp Dener     /* dL/dX += ye^T Ae */
680*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMultTransposeAdd(auglag->Ae, auglag->Ye, auglag->LgradX, auglag->LgradX));
681661095bbSAlp Dener     /* dL/dX += mu * ce^T Ae */
682*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMultTranspose(auglag->Ae, auglag->Ce, auglag->Xwork));
683*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecAXPY(auglag->LgradX, auglag->mu, auglag->Xwork));
684661095bbSAlp Dener   }
685661095bbSAlp Dener   if (tao->ineq_constrained) {
686661095bbSAlp Dener     /* compute scalar contributions */
687*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDot(auglag->Yi, auglag->Ci, &yiTcims));
688*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDot(auglag->Ci, auglag->Ci, &cimsTcims));
689661095bbSAlp Dener     /* dL/dX += yi^T Ai */
690*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMultTransposeAdd(auglag->Ai, auglag->Yi, auglag->LgradX, auglag->LgradX));
691661095bbSAlp Dener     /* dL/dX += mu * (ci - s)^T Ai */
692*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatMultTranspose(auglag->Ai, auglag->Ci, auglag->Xwork));
693*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecAXPY(auglag->LgradX, auglag->mu, auglag->Xwork));
694661095bbSAlp Dener     /* dL/dS = -[yi + mu*(ci - s)] */
695*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecWAXPY(auglag->LgradS, auglag->mu, auglag->Ci, auglag->Yi));
696*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecScale(auglag->LgradS, -1.0));
697661095bbSAlp Dener   }
698661095bbSAlp Dener   /* combine gradient together */
699*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoALMMCombinePrimal_Private(tao, auglag->LgradX, auglag->LgradS, auglag->G));
700661095bbSAlp Dener   /* compute L = f + ye^T ce + yi^T (ci - s) + 0.5*mu*||ce||^2 + 0.5*mu*||ci - s||^2 */
701661095bbSAlp Dener   auglag->Lval = auglag->fval + yeTce + yiTcims + 0.5*auglag->mu*(ceTce + cimsTcims);
702661095bbSAlp Dener   PetscFunctionReturn(0);
703661095bbSAlp Dener }
704661095bbSAlp Dener 
7057721a36fSStefano Zampini PetscErrorCode TaoALMMSubsolverObjective_Private(Tao tao, Vec P, PetscReal *Lval, void *ctx)
7067721a36fSStefano Zampini {
7077721a36fSStefano Zampini   TAO_ALMM       *auglag = (TAO_ALMM*)ctx;
7087721a36fSStefano Zampini 
7097721a36fSStefano Zampini   PetscFunctionBegin;
710*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCopy(P, auglag->P));
711*5f80ce2aSJacob Faibussowitsch   CHKERRQ((*auglag->sub_obj)(auglag->parent));
7127721a36fSStefano Zampini   *Lval = auglag->Lval;
7137721a36fSStefano Zampini   PetscFunctionReturn(0);
7147721a36fSStefano Zampini }
7157721a36fSStefano Zampini 
716661095bbSAlp Dener PetscErrorCode TaoALMMSubsolverObjectiveAndGradient_Private(Tao tao, Vec P, PetscReal *Lval, Vec G, void *ctx)
717661095bbSAlp Dener {
718661095bbSAlp Dener   TAO_ALMM       *auglag = (TAO_ALMM*)ctx;
719661095bbSAlp Dener 
720661095bbSAlp Dener   PetscFunctionBegin;
721*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCopy(P, auglag->P));
722*5f80ce2aSJacob Faibussowitsch   CHKERRQ((*auglag->sub_obj)(auglag->parent));
723*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCopy(auglag->G, G));
724661095bbSAlp Dener   *Lval = auglag->Lval;
725661095bbSAlp Dener   PetscFunctionReturn(0);
726661095bbSAlp Dener }
727