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