xref: /petsc/src/tao/pde_constrained/impls/lcl/lcl.c (revision 117ef88edefbfc12e7c19efe87a19a2e1b0acd4f)
1 #include <../src/tao/pde_constrained/impls/lcl/lcl.h>
2 static PetscErrorCode LCLComputeLagrangianAndGradient(TaoLineSearch,Vec,PetscReal*,Vec,void*);
3 static PetscErrorCode LCLComputeAugmentedLagrangianAndGradient(TaoLineSearch,Vec,PetscReal*,Vec,void*);
4 static PetscErrorCode LCLScatter(TAO_LCL*,Vec,Vec,Vec);
5 static PetscErrorCode LCLGather(TAO_LCL*,Vec,Vec,Vec);
6 
7 static PetscErrorCode TaoDestroy_LCL(Tao tao)
8 {
9   TAO_LCL        *lclP = (TAO_LCL*)tao->data;
10   PetscErrorCode ierr;
11 
12   PetscFunctionBegin;
13   if (tao->setupcalled) {
14     ierr = MatDestroy(&lclP->R);CHKERRQ(ierr);
15     ierr = VecDestroy(&lclP->lamda);CHKERRQ(ierr);
16     ierr = VecDestroy(&lclP->lamda0);CHKERRQ(ierr);
17     ierr = VecDestroy(&lclP->WL);CHKERRQ(ierr);
18     ierr = VecDestroy(&lclP->W);CHKERRQ(ierr);
19     ierr = VecDestroy(&lclP->X0);CHKERRQ(ierr);
20     ierr = VecDestroy(&lclP->G0);CHKERRQ(ierr);
21     ierr = VecDestroy(&lclP->GL);CHKERRQ(ierr);
22     ierr = VecDestroy(&lclP->GAugL);CHKERRQ(ierr);
23     ierr = VecDestroy(&lclP->dbar);CHKERRQ(ierr);
24     ierr = VecDestroy(&lclP->U);CHKERRQ(ierr);
25     ierr = VecDestroy(&lclP->U0);CHKERRQ(ierr);
26     ierr = VecDestroy(&lclP->V);CHKERRQ(ierr);
27     ierr = VecDestroy(&lclP->V0);CHKERRQ(ierr);
28     ierr = VecDestroy(&lclP->V1);CHKERRQ(ierr);
29     ierr = VecDestroy(&lclP->GU);CHKERRQ(ierr);
30     ierr = VecDestroy(&lclP->GV);CHKERRQ(ierr);
31     ierr = VecDestroy(&lclP->GU0);CHKERRQ(ierr);
32     ierr = VecDestroy(&lclP->GV0);CHKERRQ(ierr);
33     ierr = VecDestroy(&lclP->GL_U);CHKERRQ(ierr);
34     ierr = VecDestroy(&lclP->GL_V);CHKERRQ(ierr);
35     ierr = VecDestroy(&lclP->GAugL_U);CHKERRQ(ierr);
36     ierr = VecDestroy(&lclP->GAugL_V);CHKERRQ(ierr);
37     ierr = VecDestroy(&lclP->GL_U0);CHKERRQ(ierr);
38     ierr = VecDestroy(&lclP->GL_V0);CHKERRQ(ierr);
39     ierr = VecDestroy(&lclP->GAugL_U0);CHKERRQ(ierr);
40     ierr = VecDestroy(&lclP->GAugL_V0);CHKERRQ(ierr);
41     ierr = VecDestroy(&lclP->DU);CHKERRQ(ierr);
42     ierr = VecDestroy(&lclP->DV);CHKERRQ(ierr);
43     ierr = VecDestroy(&lclP->WU);CHKERRQ(ierr);
44     ierr = VecDestroy(&lclP->WV);CHKERRQ(ierr);
45     ierr = VecDestroy(&lclP->g1);CHKERRQ(ierr);
46     ierr = VecDestroy(&lclP->g2);CHKERRQ(ierr);
47     ierr = VecDestroy(&lclP->con1);CHKERRQ(ierr);
48 
49     ierr = VecDestroy(&lclP->r);CHKERRQ(ierr);
50     ierr = VecDestroy(&lclP->s);CHKERRQ(ierr);
51 
52     ierr = ISDestroy(&tao->state_is);CHKERRQ(ierr);
53     ierr = ISDestroy(&tao->design_is);CHKERRQ(ierr);
54 
55     ierr = VecScatterDestroy(&lclP->state_scatter);CHKERRQ(ierr);
56     ierr = VecScatterDestroy(&lclP->design_scatter);CHKERRQ(ierr);
57   }
58   ierr = MatDestroy(&lclP->R);CHKERRQ(ierr);
59   ierr = PetscFree(tao->data);CHKERRQ(ierr);
60   PetscFunctionReturn(0);
61 }
62 
63 static PetscErrorCode TaoSetFromOptions_LCL(PetscOptionItems *PetscOptionsObject,Tao tao)
64 {
65   TAO_LCL        *lclP = (TAO_LCL*)tao->data;
66   PetscErrorCode ierr;
67 
68   PetscFunctionBegin;
69   ierr = PetscOptionsHead(PetscOptionsObject,"Linearly-Constrained Augmented Lagrangian Method for PDE-constrained optimization");CHKERRQ(ierr);
70   ierr = PetscOptionsReal("-tao_lcl_eps1","epsilon 1 tolerance","",lclP->eps1,&lclP->eps1,NULL);CHKERRQ(ierr);
71   ierr = PetscOptionsReal("-tao_lcl_eps2","epsilon 2 tolerance","",lclP->eps2,&lclP->eps2,NULL);CHKERRQ(ierr);
72   ierr = PetscOptionsReal("-tao_lcl_rho0","init value for rho","",lclP->rho0,&lclP->rho0,NULL);CHKERRQ(ierr);
73   ierr = PetscOptionsReal("-tao_lcl_rhomax","max value for rho","",lclP->rhomax,&lclP->rhomax,NULL);CHKERRQ(ierr);
74   lclP->phase2_niter = 1;
75   ierr = PetscOptionsInt("-tao_lcl_phase2_niter","Number of phase 2 iterations in LCL algorithm","",lclP->phase2_niter,&lclP->phase2_niter,NULL);CHKERRQ(ierr);
76   lclP->verbose = PETSC_FALSE;
77   ierr = PetscOptionsBool("-tao_lcl_verbose","Print verbose output","",lclP->verbose,&lclP->verbose,NULL);CHKERRQ(ierr);
78   lclP->tau[0] = lclP->tau[1] = lclP->tau[2] = lclP->tau[3] = 1.0e-4;
79   ierr = PetscOptionsReal("-tao_lcl_tola","Tolerance for first forward solve","",lclP->tau[0],&lclP->tau[0],NULL);CHKERRQ(ierr);
80   ierr = PetscOptionsReal("-tao_lcl_tolb","Tolerance for first adjoint solve","",lclP->tau[1],&lclP->tau[1],NULL);CHKERRQ(ierr);
81   ierr = PetscOptionsReal("-tao_lcl_tolc","Tolerance for second forward solve","",lclP->tau[2],&lclP->tau[2],NULL);CHKERRQ(ierr);
82   ierr = PetscOptionsReal("-tao_lcl_told","Tolerance for second adjoint solve","",lclP->tau[3],&lclP->tau[3],NULL);CHKERRQ(ierr);
83   ierr = PetscOptionsTail();CHKERRQ(ierr);
84   ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
85   ierr = MatSetFromOptions(lclP->R);CHKERRQ(ierr);
86   PetscFunctionReturn(0);
87 }
88 
89 static PetscErrorCode TaoView_LCL(Tao tao, PetscViewer viewer)
90 {
91   return 0;
92 }
93 
94 static PetscErrorCode TaoSetup_LCL(Tao tao)
95 {
96   TAO_LCL        *lclP = (TAO_LCL*)tao->data;
97   PetscInt       lo, hi, nlocalstate, nlocaldesign;
98   PetscErrorCode ierr;
99   IS             is_state, is_design;
100 
101   PetscFunctionBegin;
102   if (!tao->state_is) SETERRQ(PetscObjectComm((PetscObject)tao),PETSC_ERR_ARG_WRONGSTATE,"LCL Solver requires an initial state index set -- use TaoSetStateIS()");
103   ierr = VecDuplicate(tao->solution, &tao->gradient);CHKERRQ(ierr);
104   ierr = VecDuplicate(tao->solution, &tao->stepdirection);CHKERRQ(ierr);
105   ierr = VecDuplicate(tao->solution, &lclP->W);CHKERRQ(ierr);
106   ierr = VecDuplicate(tao->solution, &lclP->X0);CHKERRQ(ierr);
107   ierr = VecDuplicate(tao->solution, &lclP->G0);CHKERRQ(ierr);
108   ierr = VecDuplicate(tao->solution, &lclP->GL);CHKERRQ(ierr);
109   ierr = VecDuplicate(tao->solution, &lclP->GAugL);CHKERRQ(ierr);
110 
111   ierr = VecDuplicate(tao->constraints, &lclP->lamda);CHKERRQ(ierr);
112   ierr = VecDuplicate(tao->constraints, &lclP->WL);CHKERRQ(ierr);
113   ierr = VecDuplicate(tao->constraints, &lclP->lamda0);CHKERRQ(ierr);
114   ierr = VecDuplicate(tao->constraints, &lclP->con1);CHKERRQ(ierr);
115 
116   ierr = VecSet(lclP->lamda,0.0);CHKERRQ(ierr);
117 
118   ierr = VecGetSize(tao->solution, &lclP->n);CHKERRQ(ierr);
119   ierr = VecGetSize(tao->constraints, &lclP->m);CHKERRQ(ierr);
120 
121   ierr = VecCreate(((PetscObject)tao)->comm,&lclP->U);CHKERRQ(ierr);
122   ierr = VecCreate(((PetscObject)tao)->comm,&lclP->V);CHKERRQ(ierr);
123   ierr = ISGetLocalSize(tao->state_is,&nlocalstate);CHKERRQ(ierr);
124   ierr = ISGetLocalSize(tao->design_is,&nlocaldesign);CHKERRQ(ierr);
125   ierr = VecSetSizes(lclP->U,nlocalstate,lclP->m);CHKERRQ(ierr);
126   ierr = VecSetSizes(lclP->V,nlocaldesign,lclP->n-lclP->m);CHKERRQ(ierr);
127   ierr = VecSetType(lclP->U,((PetscObject)(tao->solution))->type_name);CHKERRQ(ierr);
128   ierr = VecSetType(lclP->V,((PetscObject)(tao->solution))->type_name);CHKERRQ(ierr);
129   ierr = VecSetFromOptions(lclP->U);CHKERRQ(ierr);
130   ierr = VecSetFromOptions(lclP->V);CHKERRQ(ierr);
131   ierr = VecDuplicate(lclP->U,&lclP->DU);CHKERRQ(ierr);
132   ierr = VecDuplicate(lclP->U,&lclP->U0);CHKERRQ(ierr);
133   ierr = VecDuplicate(lclP->U,&lclP->GU);CHKERRQ(ierr);
134   ierr = VecDuplicate(lclP->U,&lclP->GU0);CHKERRQ(ierr);
135   ierr = VecDuplicate(lclP->U,&lclP->GAugL_U);CHKERRQ(ierr);
136   ierr = VecDuplicate(lclP->U,&lclP->GL_U);CHKERRQ(ierr);
137   ierr = VecDuplicate(lclP->U,&lclP->GAugL_U0);CHKERRQ(ierr);
138   ierr = VecDuplicate(lclP->U,&lclP->GL_U0);CHKERRQ(ierr);
139   ierr = VecDuplicate(lclP->U,&lclP->WU);CHKERRQ(ierr);
140   ierr = VecDuplicate(lclP->U,&lclP->r);CHKERRQ(ierr);
141   ierr = VecDuplicate(lclP->V,&lclP->V0);CHKERRQ(ierr);
142   ierr = VecDuplicate(lclP->V,&lclP->V1);CHKERRQ(ierr);
143   ierr = VecDuplicate(lclP->V,&lclP->DV);CHKERRQ(ierr);
144   ierr = VecDuplicate(lclP->V,&lclP->s);CHKERRQ(ierr);
145   ierr = VecDuplicate(lclP->V,&lclP->GV);CHKERRQ(ierr);
146   ierr = VecDuplicate(lclP->V,&lclP->GV0);CHKERRQ(ierr);
147   ierr = VecDuplicate(lclP->V,&lclP->dbar);CHKERRQ(ierr);
148   ierr = VecDuplicate(lclP->V,&lclP->GAugL_V);CHKERRQ(ierr);
149   ierr = VecDuplicate(lclP->V,&lclP->GL_V);CHKERRQ(ierr);
150   ierr = VecDuplicate(lclP->V,&lclP->GAugL_V0);CHKERRQ(ierr);
151   ierr = VecDuplicate(lclP->V,&lclP->GL_V0);CHKERRQ(ierr);
152   ierr = VecDuplicate(lclP->V,&lclP->WV);CHKERRQ(ierr);
153   ierr = VecDuplicate(lclP->V,&lclP->g1);CHKERRQ(ierr);
154   ierr = VecDuplicate(lclP->V,&lclP->g2);CHKERRQ(ierr);
155 
156   /* create scatters for state, design subvecs */
157   ierr = VecGetOwnershipRange(lclP->U,&lo,&hi);CHKERRQ(ierr);
158   ierr = ISCreateStride(((PetscObject)lclP->U)->comm,hi-lo,lo,1,&is_state);CHKERRQ(ierr);
159   ierr = VecGetOwnershipRange(lclP->V,&lo,&hi);CHKERRQ(ierr);
160   if (0) {
161     PetscInt sizeU,sizeV;
162     ierr = VecGetSize(lclP->U,&sizeU);CHKERRQ(ierr);
163     ierr = VecGetSize(lclP->V,&sizeV);CHKERRQ(ierr);
164     ierr = PetscPrintf(PETSC_COMM_WORLD,"size(U)=%D, size(V)=%D\n",sizeU,sizeV);CHKERRQ(ierr);
165   }
166   ierr = ISCreateStride(((PetscObject)lclP->V)->comm,hi-lo,lo,1,&is_design);CHKERRQ(ierr);
167   ierr = VecScatterCreate(tao->solution,tao->state_is,lclP->U,is_state,&lclP->state_scatter);CHKERRQ(ierr);
168   ierr = VecScatterCreate(tao->solution,tao->design_is,lclP->V,is_design,&lclP->design_scatter);CHKERRQ(ierr);
169   ierr = ISDestroy(&is_state);CHKERRQ(ierr);
170   ierr = ISDestroy(&is_design);CHKERRQ(ierr);
171   PetscFunctionReturn(0);
172 }
173 
174 static PetscErrorCode TaoSolve_LCL(Tao tao)
175 {
176   TAO_LCL                      *lclP = (TAO_LCL*)tao->data;
177   PetscInt                     phase2_iter,nlocal,its;
178   TaoLineSearchConvergedReason ls_reason = TAOLINESEARCH_CONTINUE_ITERATING;
179   PetscReal                    step=1.0,f, descent, aldescent;
180   PetscReal                    cnorm, mnorm;
181   PetscReal                    adec,r2,rGL_U,rWU;
182   PetscBool                    set,pset,flag,pflag,symmetric;
183   PetscErrorCode               ierr;
184 
185   PetscFunctionBegin;
186   lclP->rho = lclP->rho0;
187   ierr = VecGetLocalSize(lclP->U,&nlocal);CHKERRQ(ierr);
188   ierr = VecGetLocalSize(lclP->V,&nlocal);CHKERRQ(ierr);
189   ierr = MatSetSizes(lclP->R, nlocal, nlocal, lclP->n-lclP->m, lclP->n-lclP->m);CHKERRQ(ierr);
190   ierr = MatLMVMAllocate(lclP->R,lclP->V,lclP->V);CHKERRQ(ierr);
191   lclP->recompute_jacobian_flag = PETSC_TRUE;
192 
193   /* Scatter to U,V */
194   ierr = LCLScatter(lclP,tao->solution,lclP->U,lclP->V);CHKERRQ(ierr);
195 
196   /* Evaluate Function, Gradient, Constraints, and Jacobian */
197   ierr = TaoComputeObjectiveAndGradient(tao,tao->solution,&f,tao->gradient);CHKERRQ(ierr);
198   ierr = TaoComputeJacobianState(tao,tao->solution,tao->jacobian_state,tao->jacobian_state_pre,tao->jacobian_state_inv);CHKERRQ(ierr);
199   ierr = TaoComputeJacobianDesign(tao,tao->solution,tao->jacobian_design);CHKERRQ(ierr);
200   ierr = TaoComputeConstraints(tao,tao->solution, tao->constraints);CHKERRQ(ierr);
201 
202   /* Scatter gradient to GU,GV */
203   ierr = LCLScatter(lclP,tao->gradient,lclP->GU,lclP->GV);CHKERRQ(ierr);
204 
205   /* Evaluate Lagrangian function and gradient */
206   /* p0 */
207   ierr = VecSet(lclP->lamda,0.0);CHKERRQ(ierr); /*  Initial guess in CG */
208   ierr = MatIsSymmetricKnown(tao->jacobian_state,&set,&flag);CHKERRQ(ierr);
209   if (tao->jacobian_state_pre) {
210     ierr = MatIsSymmetricKnown(tao->jacobian_state_pre,&pset,&pflag);CHKERRQ(ierr);
211   } else {
212     pset = pflag = PETSC_TRUE;
213   }
214   if (set && pset && flag && pflag) symmetric = PETSC_TRUE;
215   else symmetric = PETSC_FALSE;
216 
217   lclP->solve_type = LCL_ADJOINT2;
218   if (tao->jacobian_state_inv) {
219     if (symmetric) {
220       ierr = MatMult(tao->jacobian_state_inv, lclP->GU, lclP->lamda);CHKERRQ(ierr); } else {
221       ierr = MatMultTranspose(tao->jacobian_state_inv, lclP->GU, lclP->lamda);CHKERRQ(ierr);
222     }
223   } else {
224     ierr = KSPSetOperators(tao->ksp, tao->jacobian_state, tao->jacobian_state_pre);CHKERRQ(ierr);
225     if (symmetric) {
226       ierr = KSPSolve(tao->ksp, lclP->GU,  lclP->lamda);CHKERRQ(ierr);
227     } else {
228       ierr = KSPSolveTranspose(tao->ksp, lclP->GU,  lclP->lamda);CHKERRQ(ierr);
229     }
230     ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr);
231     tao->ksp_its+=its;
232     tao->ksp_tot_its+=its;
233   }
234   ierr = VecCopy(lclP->lamda,lclP->lamda0);CHKERRQ(ierr);
235   ierr = LCLComputeAugmentedLagrangianAndGradient(tao->linesearch,tao->solution,&lclP->aug,lclP->GAugL,tao);CHKERRQ(ierr);
236 
237   ierr = LCLScatter(lclP,lclP->GL,lclP->GL_U,lclP->GL_V);CHKERRQ(ierr);
238   ierr = LCLScatter(lclP,lclP->GAugL,lclP->GAugL_U,lclP->GAugL_V);CHKERRQ(ierr);
239 
240   /* Evaluate constraint norm */
241   ierr = VecNorm(tao->constraints, NORM_2, &cnorm);CHKERRQ(ierr);
242   ierr = VecNorm(lclP->GAugL, NORM_2, &mnorm);CHKERRQ(ierr);
243 
244   /* Monitor convergence */
245   tao->reason = TAO_CONTINUE_ITERATING;
246   ierr = TaoLogConvergenceHistory(tao,f,mnorm,cnorm,tao->ksp_its);CHKERRQ(ierr);
247   ierr = TaoMonitor(tao,tao->niter,f,mnorm,cnorm,step);CHKERRQ(ierr);
248   ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
249 
250   while (tao->reason == TAO_CONTINUE_ITERATING) {
251     /* Call general purpose update function */
252     if (tao->ops->update) {
253       ierr = (*tao->ops->update)(tao, tao->niter, tao->user_update);CHKERRQ(ierr);
254     }
255     tao->ksp_its=0;
256     /* Compute a descent direction for the linearly constrained subproblem
257        minimize f(u+du, v+dv)
258        s.t. A(u0,v0)du + B(u0,v0)dv = -g(u0,v0) */
259 
260     /* Store the points around the linearization */
261     ierr = VecCopy(lclP->U, lclP->U0);CHKERRQ(ierr);
262     ierr = VecCopy(lclP->V, lclP->V0);CHKERRQ(ierr);
263     ierr = VecCopy(lclP->GU,lclP->GU0);CHKERRQ(ierr);
264     ierr = VecCopy(lclP->GV,lclP->GV0);CHKERRQ(ierr);
265     ierr = VecCopy(lclP->GAugL_U,lclP->GAugL_U0);CHKERRQ(ierr);
266     ierr = VecCopy(lclP->GAugL_V,lclP->GAugL_V0);CHKERRQ(ierr);
267     ierr = VecCopy(lclP->GL_U,lclP->GL_U0);CHKERRQ(ierr);
268     ierr = VecCopy(lclP->GL_V,lclP->GL_V0);CHKERRQ(ierr);
269 
270     lclP->aug0 = lclP->aug;
271     lclP->lgn0 = lclP->lgn;
272 
273     /* Given the design variables, we need to project the current iterate
274        onto the linearized constraint.  We choose to fix the design variables
275        and solve the linear system for the state variables.  The resulting
276        point is the Newton direction */
277 
278     /* Solve r = A\con */
279     lclP->solve_type = LCL_FORWARD1;
280     ierr = VecSet(lclP->r,0.0);CHKERRQ(ierr); /*  Initial guess in CG */
281 
282     if (tao->jacobian_state_inv) {
283       ierr = MatMult(tao->jacobian_state_inv, tao->constraints, lclP->r);CHKERRQ(ierr);
284     } else {
285       ierr = KSPSetOperators(tao->ksp, tao->jacobian_state, tao->jacobian_state_pre);CHKERRQ(ierr);
286       ierr = KSPSolve(tao->ksp, tao->constraints,  lclP->r);CHKERRQ(ierr);
287       ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr);
288       tao->ksp_its+=its;
289       tao->ksp_tot_its+=tao->ksp_its;
290     }
291 
292     /* Set design step direction dv to zero */
293     ierr = VecSet(lclP->s, 0.0);CHKERRQ(ierr);
294 
295     /*
296        Check sufficient descent for constraint merit function .5*||con||^2
297        con' Ak r >= eps1 ||r||^(2+eps2)
298     */
299 
300     /* Compute WU= Ak' * con */
301     if (symmetric)  {
302       ierr = MatMult(tao->jacobian_state,tao->constraints,lclP->WU);CHKERRQ(ierr);
303     } else {
304       ierr = MatMultTranspose(tao->jacobian_state,tao->constraints,lclP->WU);CHKERRQ(ierr);
305     }
306     /* Compute r * Ak' * con */
307     ierr = VecDot(lclP->r,lclP->WU,&rWU);CHKERRQ(ierr);
308 
309     /* compute ||r||^(2+eps2) */
310     ierr = VecNorm(lclP->r,NORM_2,&r2);CHKERRQ(ierr);
311     r2 = PetscPowScalar(r2,2.0+lclP->eps2);
312     adec = lclP->eps1 * r2;
313 
314     if (rWU < adec) {
315       ierr = PetscInfo(tao,"Newton direction not descent for constraint, feasibility phase required\n");CHKERRQ(ierr);
316       if (lclP->verbose) {
317         ierr = PetscPrintf(PETSC_COMM_WORLD,"Newton direction not descent for constraint: %g -- using steepest descent\n",(double)descent);CHKERRQ(ierr);
318       }
319 
320       ierr = PetscInfo(tao,"Using steepest descent direction instead.\n");CHKERRQ(ierr);
321       ierr = VecSet(lclP->r,0.0);CHKERRQ(ierr);
322       ierr = VecAXPY(lclP->r,-1.0,lclP->WU);CHKERRQ(ierr);
323       ierr = VecDot(lclP->r,lclP->r,&rWU);CHKERRQ(ierr);
324       ierr = VecNorm(lclP->r,NORM_2,&r2);CHKERRQ(ierr);
325       r2 = PetscPowScalar(r2,2.0+lclP->eps2);
326       ierr = VecDot(lclP->r,lclP->GAugL_U,&descent);CHKERRQ(ierr);
327       adec = lclP->eps1 * r2;
328     }
329 
330 
331     /*
332        Check descent for aug. lagrangian
333        r' (GUk - Ak'*yk - rho*Ak'*con) <= -eps1 ||r||^(2+eps2)
334           GL_U = GUk - Ak'*yk
335           WU   = Ak'*con
336                                          adec=eps1||r||^(2+eps2)
337 
338        ==>
339        Check r'GL_U - rho*r'WU <= adec
340     */
341 
342     ierr = VecDot(lclP->r,lclP->GL_U,&rGL_U);CHKERRQ(ierr);
343     aldescent =  rGL_U - lclP->rho*rWU;
344     if (aldescent > -adec) {
345       if (lclP->verbose) {
346         ierr = PetscPrintf(PETSC_COMM_WORLD," Newton direction not descent for augmented Lagrangian: %g",(double)aldescent);CHKERRQ(ierr);
347       }
348       ierr = PetscInfo1(tao,"Newton direction not descent for augmented Lagrangian: %g\n",(double)aldescent);CHKERRQ(ierr);
349       lclP->rho =  (rGL_U - adec)/rWU;
350       if (lclP->rho > lclP->rhomax) {
351         lclP->rho = lclP->rhomax;
352         SETERRQ1(PETSC_COMM_WORLD,0,"rho=%g > rhomax, case not implemented.  Increase rhomax (-tao_lcl_rhomax)",(double)lclP->rho);
353       }
354       if (lclP->verbose) {
355         ierr = PetscPrintf(PETSC_COMM_WORLD,"  Increasing penalty parameter to %g\n",(double)lclP->rho);CHKERRQ(ierr);
356       }
357       ierr = PetscInfo1(tao,"  Increasing penalty parameter to %g\n",(double)lclP->rho);CHKERRQ(ierr);
358     }
359 
360     ierr = LCLComputeAugmentedLagrangianAndGradient(tao->linesearch,tao->solution,&lclP->aug,lclP->GAugL,tao);CHKERRQ(ierr);
361     ierr = LCLScatter(lclP,lclP->GL,lclP->GL_U,lclP->GL_V);CHKERRQ(ierr);
362     ierr = LCLScatter(lclP,lclP->GAugL,lclP->GAugL_U,lclP->GAugL_V);CHKERRQ(ierr);
363 
364     /* We now minimize the augmented Lagrangian along the Newton direction */
365     ierr = VecScale(lclP->r,-1.0);CHKERRQ(ierr);
366     ierr = LCLGather(lclP, lclP->r,lclP->s,tao->stepdirection);CHKERRQ(ierr);
367     ierr = VecScale(lclP->r,-1.0);CHKERRQ(ierr);
368     ierr = LCLGather(lclP, lclP->GAugL_U0, lclP->GAugL_V0, lclP->GAugL);CHKERRQ(ierr);
369     ierr = LCLGather(lclP, lclP->U0,lclP->V0,lclP->X0);CHKERRQ(ierr);
370 
371     lclP->recompute_jacobian_flag = PETSC_TRUE;
372 
373     ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);CHKERRQ(ierr);
374     ierr = TaoLineSearchSetType(tao->linesearch, TAOLINESEARCHMT);CHKERRQ(ierr);
375     ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
376     ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &lclP->aug, lclP->GAugL, tao->stepdirection, &step, &ls_reason);CHKERRQ(ierr);
377     if (lclP->verbose) {
378       ierr = PetscPrintf(PETSC_COMM_WORLD,"Steplength = %g\n",(double)step);CHKERRQ(ierr);
379     }
380 
381     ierr = LCLScatter(lclP,tao->solution,lclP->U,lclP->V);CHKERRQ(ierr);
382     ierr = TaoComputeObjectiveAndGradient(tao,tao->solution,&f,tao->gradient);CHKERRQ(ierr);
383     ierr = LCLScatter(lclP,tao->gradient,lclP->GU,lclP->GV);CHKERRQ(ierr);
384 
385     ierr = LCLScatter(lclP,lclP->GAugL,lclP->GAugL_U,lclP->GAugL_V);CHKERRQ(ierr);
386 
387     /* Check convergence */
388     ierr = VecNorm(lclP->GAugL, NORM_2, &mnorm);CHKERRQ(ierr);
389     ierr = VecNorm(tao->constraints, NORM_2, &cnorm);CHKERRQ(ierr);
390     ierr = TaoLogConvergenceHistory(tao,f,mnorm,cnorm,tao->ksp_its);CHKERRQ(ierr);
391     ierr = TaoMonitor(tao,tao->niter,f,mnorm,cnorm,step);CHKERRQ(ierr);
392     ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
393     if (tao->reason != TAO_CONTINUE_ITERATING){
394       break;
395     }
396 
397     /* TODO: use a heuristic to choose how many iterations should be performed within phase 2 */
398     for (phase2_iter=0; phase2_iter<lclP->phase2_niter; phase2_iter++){
399       /* We now minimize the objective function starting from the fraction of
400          the Newton point accepted by applying one step of a reduced-space
401          method.  The optimization problem is:
402 
403          minimize f(u+du, v+dv)
404          s. t.    A(u0,v0)du + B(u0,v0)du = -alpha g(u0,v0)
405 
406          In particular, we have that
407          du = -inv(A)*(Bdv + alpha g) */
408 
409       ierr = TaoComputeJacobianState(tao,lclP->X0,tao->jacobian_state,tao->jacobian_state_pre,tao->jacobian_state_inv);CHKERRQ(ierr);
410       ierr = TaoComputeJacobianDesign(tao,lclP->X0,tao->jacobian_design);CHKERRQ(ierr);
411 
412       /* Store V and constraints */
413       ierr = VecCopy(lclP->V, lclP->V1);CHKERRQ(ierr);
414       ierr = VecCopy(tao->constraints,lclP->con1);CHKERRQ(ierr);
415 
416       /* Compute multipliers */
417       /* p1 */
418       ierr = VecSet(lclP->lamda,0.0);CHKERRQ(ierr); /*  Initial guess in CG */
419       lclP->solve_type = LCL_ADJOINT1;
420       ierr = MatIsSymmetricKnown(tao->jacobian_state,&set,&flag);CHKERRQ(ierr);
421       if (tao->jacobian_state_pre) {
422         ierr = MatIsSymmetricKnown(tao->jacobian_state_pre,&pset,&pflag);CHKERRQ(ierr);
423       } else {
424         pset = pflag = PETSC_TRUE;
425       }
426       if (set && pset && flag && pflag) symmetric = PETSC_TRUE;
427       else symmetric = PETSC_FALSE;
428 
429       if (tao->jacobian_state_inv) {
430         if (symmetric) {
431           ierr = MatMult(tao->jacobian_state_inv, lclP->GAugL_U, lclP->lamda);CHKERRQ(ierr);
432         } else {
433           ierr = MatMultTranspose(tao->jacobian_state_inv, lclP->GAugL_U, lclP->lamda);CHKERRQ(ierr);
434         }
435       } else {
436         if (symmetric) {
437           ierr = KSPSolve(tao->ksp, lclP->GAugL_U,  lclP->lamda);CHKERRQ(ierr);
438         } else {
439           ierr = KSPSolveTranspose(tao->ksp, lclP->GAugL_U,  lclP->lamda);CHKERRQ(ierr);
440         }
441         ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr);
442         tao->ksp_its+=its;
443         tao->ksp_tot_its+=its;
444       }
445       ierr = MatMultTranspose(tao->jacobian_design,lclP->lamda,lclP->g1);CHKERRQ(ierr);
446       ierr = VecAXPY(lclP->g1,-1.0,lclP->GAugL_V);CHKERRQ(ierr);
447 
448       /* Compute the limited-memory quasi-newton direction */
449       if (tao->niter > 0) {
450         ierr = MatSolve(lclP->R,lclP->g1,lclP->s);CHKERRQ(ierr);
451         ierr = VecDot(lclP->s,lclP->g1,&descent);CHKERRQ(ierr);
452         if (descent <= 0) {
453           if (lclP->verbose) {
454             ierr = PetscPrintf(PETSC_COMM_WORLD,"Reduced-space direction not descent: %g\n",(double)descent);CHKERRQ(ierr);
455           }
456           ierr = VecCopy(lclP->g1,lclP->s);CHKERRQ(ierr);
457         }
458       } else {
459         ierr = VecCopy(lclP->g1,lclP->s);CHKERRQ(ierr);
460       }
461       ierr = VecScale(lclP->g1,-1.0);CHKERRQ(ierr);
462 
463       /* Recover the full space direction */
464       ierr = MatMult(tao->jacobian_design,lclP->s,lclP->WU);CHKERRQ(ierr);
465       /* ierr = VecSet(lclP->r,0.0);CHKERRQ(ierr); */ /*  Initial guess in CG */
466       lclP->solve_type = LCL_FORWARD2;
467       if (tao->jacobian_state_inv) {
468         ierr = MatMult(tao->jacobian_state_inv,lclP->WU,lclP->r);CHKERRQ(ierr);
469       } else {
470         ierr = KSPSolve(tao->ksp, lclP->WU, lclP->r);CHKERRQ(ierr);
471         ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr);
472         tao->ksp_its += its;
473         tao->ksp_tot_its+=its;
474       }
475 
476       /* We now minimize the augmented Lagrangian along the direction -r,s */
477       ierr = VecScale(lclP->r, -1.0);CHKERRQ(ierr);
478       ierr = LCLGather(lclP,lclP->r,lclP->s,tao->stepdirection);CHKERRQ(ierr);
479       ierr = VecScale(lclP->r, -1.0);CHKERRQ(ierr);
480       lclP->recompute_jacobian_flag = PETSC_TRUE;
481 
482       ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);CHKERRQ(ierr);
483       ierr = TaoLineSearchSetType(tao->linesearch,TAOLINESEARCHMT);CHKERRQ(ierr);
484       ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
485       ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &lclP->aug, lclP->GAugL, tao->stepdirection,&step,&ls_reason);CHKERRQ(ierr);
486       if (lclP->verbose){
487         ierr = PetscPrintf(PETSC_COMM_WORLD,"Reduced-space steplength =  %g\n",(double)step);CHKERRQ(ierr);
488       }
489 
490       ierr = LCLScatter(lclP,tao->solution,lclP->U,lclP->V);CHKERRQ(ierr);
491       ierr = LCLScatter(lclP,lclP->GL,lclP->GL_U,lclP->GL_V);CHKERRQ(ierr);
492       ierr = LCLScatter(lclP,lclP->GAugL,lclP->GAugL_U,lclP->GAugL_V);CHKERRQ(ierr);
493       ierr = TaoComputeObjectiveAndGradient(tao,tao->solution,&f,tao->gradient);CHKERRQ(ierr);
494       ierr = LCLScatter(lclP,tao->gradient,lclP->GU,lclP->GV);CHKERRQ(ierr);
495 
496       /* Compute the reduced gradient at the new point */
497 
498       ierr = TaoComputeJacobianState(tao,lclP->X0,tao->jacobian_state,tao->jacobian_state_pre,tao->jacobian_state_inv);CHKERRQ(ierr);
499       ierr = TaoComputeJacobianDesign(tao,lclP->X0,tao->jacobian_design);CHKERRQ(ierr);
500 
501       /* p2 */
502       /* Compute multipliers, using lamda-rho*con as an initial guess in PCG */
503       if (phase2_iter==0){
504         ierr = VecWAXPY(lclP->lamda,-lclP->rho,lclP->con1,lclP->lamda0);CHKERRQ(ierr);
505       } else {
506         ierr = VecAXPY(lclP->lamda,-lclP->rho,tao->constraints);CHKERRQ(ierr);
507       }
508 
509       ierr = MatIsSymmetricKnown(tao->jacobian_state,&set,&flag);CHKERRQ(ierr);
510       if (tao->jacobian_state_pre) {
511         ierr = MatIsSymmetricKnown(tao->jacobian_state_pre,&pset,&pflag);CHKERRQ(ierr);
512       } else {
513         pset = pflag = PETSC_TRUE;
514       }
515       if (set && pset && flag && pflag) symmetric = PETSC_TRUE;
516       else symmetric = PETSC_FALSE;
517 
518       lclP->solve_type = LCL_ADJOINT2;
519       if (tao->jacobian_state_inv) {
520         if (symmetric) {
521           ierr = MatMult(tao->jacobian_state_inv, lclP->GU, lclP->lamda);CHKERRQ(ierr);
522         } else {
523           ierr = MatMultTranspose(tao->jacobian_state_inv, lclP->GU, lclP->lamda);CHKERRQ(ierr);
524         }
525       } else {
526         if (symmetric) {
527           ierr = KSPSolve(tao->ksp, lclP->GU,  lclP->lamda);CHKERRQ(ierr);
528         } else {
529           ierr = KSPSolveTranspose(tao->ksp, lclP->GU,  lclP->lamda);CHKERRQ(ierr);
530         }
531         ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr);
532         tao->ksp_its += its;
533         tao->ksp_tot_its += its;
534       }
535 
536       ierr = MatMultTranspose(tao->jacobian_design,lclP->lamda,lclP->g2);CHKERRQ(ierr);
537       ierr = VecAXPY(lclP->g2,-1.0,lclP->GV);CHKERRQ(ierr);
538 
539       ierr = VecScale(lclP->g2,-1.0);CHKERRQ(ierr);
540 
541       /* Update the quasi-newton approximation */
542       ierr = MatLMVMUpdate(lclP->R,lclP->V,lclP->g2);CHKERRQ(ierr);
543       /* Use "-tao_ls_type gpcg -tao_ls_ftol 0 -tao_lmm_broyden_phi 0.0 -tao_lmm_scale_type scalar" to obtain agreement with Matlab code */
544 
545     }
546 
547     ierr = VecCopy(lclP->lamda,lclP->lamda0);CHKERRQ(ierr);
548 
549     /* Evaluate Function, Gradient, Constraints, and Jacobian */
550     ierr = TaoComputeObjectiveAndGradient(tao,tao->solution,&f,tao->gradient);CHKERRQ(ierr);
551     ierr = LCLScatter(lclP,tao->solution,lclP->U,lclP->V);CHKERRQ(ierr);
552     ierr = LCLScatter(lclP,tao->gradient,lclP->GU,lclP->GV);CHKERRQ(ierr);
553 
554     ierr = TaoComputeJacobianState(tao,tao->solution,tao->jacobian_state,tao->jacobian_state_pre,tao->jacobian_state_inv);CHKERRQ(ierr);
555     ierr = TaoComputeJacobianDesign(tao,tao->solution,tao->jacobian_design);CHKERRQ(ierr);
556     ierr = TaoComputeConstraints(tao,tao->solution, tao->constraints);CHKERRQ(ierr);
557 
558     ierr = LCLComputeAugmentedLagrangianAndGradient(tao->linesearch,tao->solution,&lclP->aug,lclP->GAugL,tao);CHKERRQ(ierr);
559 
560     ierr = VecNorm(lclP->GAugL, NORM_2, &mnorm);CHKERRQ(ierr);
561 
562     /* Evaluate constraint norm */
563     ierr = VecNorm(tao->constraints, NORM_2, &cnorm);CHKERRQ(ierr);
564 
565     /* Monitor convergence */
566     tao->niter++;
567     ierr = TaoLogConvergenceHistory(tao,f,mnorm,cnorm,tao->ksp_its);CHKERRQ(ierr);
568     ierr = TaoMonitor(tao,tao->niter,f,mnorm,cnorm,step);CHKERRQ(ierr);
569     ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
570   }
571   PetscFunctionReturn(0);
572 }
573 
574 /*MC
575  TAOLCL - linearly constrained lagrangian method for pde-constrained optimization
576 
577 + -tao_lcl_eps1 - epsilon 1 tolerance
578 . -tao_lcl_eps2","epsilon 2 tolerance","",lclP->eps2,&lclP->eps2,NULL);CHKERRQ(ierr);
579 . -tao_lcl_rho0","init value for rho","",lclP->rho0,&lclP->rho0,NULL);CHKERRQ(ierr);
580 . -tao_lcl_rhomax","max value for rho","",lclP->rhomax,&lclP->rhomax,NULL);CHKERRQ(ierr);
581 . -tao_lcl_phase2_niter - Number of phase 2 iterations in LCL algorithm
582 . -tao_lcl_verbose - Print verbose output if True
583 . -tao_lcl_tola - Tolerance for first forward solve
584 . -tao_lcl_tolb - Tolerance for first adjoint solve
585 . -tao_lcl_tolc - Tolerance for second forward solve
586 - -tao_lcl_told - Tolerance for second adjoint solve
587 
588   Level: beginner
589 M*/
590 PETSC_EXTERN PetscErrorCode TaoCreate_LCL(Tao tao)
591 {
592   TAO_LCL        *lclP;
593   PetscErrorCode ierr;
594   const char     *morethuente_type = TAOLINESEARCHMT;
595 
596   PetscFunctionBegin;
597   tao->ops->setup = TaoSetup_LCL;
598   tao->ops->solve = TaoSolve_LCL;
599   tao->ops->view = TaoView_LCL;
600   tao->ops->setfromoptions = TaoSetFromOptions_LCL;
601   tao->ops->destroy = TaoDestroy_LCL;
602   ierr = PetscNewLog(tao,&lclP);CHKERRQ(ierr);
603   tao->data = (void*)lclP;
604 
605   /* Override default settings (unless already changed) */
606   if (!tao->max_it_changed) tao->max_it = 200;
607   if (!tao->catol_changed) tao->catol = 1.0e-4;
608   if (!tao->gatol_changed) tao->gttol = 1.0e-4;
609   if (!tao->grtol_changed) tao->gttol = 1.0e-4;
610   if (!tao->gttol_changed) tao->gttol = 1.0e-4;
611   lclP->rho0 = 1.0e-4;
612   lclP->rhomax=1e5;
613   lclP->eps1 = 1.0e-8;
614   lclP->eps2 = 0.0;
615   lclP->solve_type=2;
616   lclP->tau[0] = lclP->tau[1] = lclP->tau[2] = lclP->tau[3] = 1.0e-4;
617   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);CHKERRQ(ierr);
618   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr);
619   ierr = TaoLineSearchSetType(tao->linesearch, morethuente_type);CHKERRQ(ierr);
620   ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr);
621 
622   ierr = TaoLineSearchSetObjectiveAndGradientRoutine(tao->linesearch,LCLComputeAugmentedLagrangianAndGradient, tao);CHKERRQ(ierr);
623   ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr);
624   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1);CHKERRQ(ierr);
625   ierr = KSPSetOptionsPrefix(tao->ksp, tao->hdr.prefix);CHKERRQ(ierr);
626   ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
627 
628   ierr = MatCreate(((PetscObject)tao)->comm, &lclP->R);CHKERRQ(ierr);
629   ierr = MatSetType(lclP->R, MATLMVMBFGS);CHKERRQ(ierr);
630   PetscFunctionReturn(0);
631 }
632 
633 static PetscErrorCode LCLComputeLagrangianAndGradient(TaoLineSearch ls, Vec X, PetscReal *f, Vec G, void *ptr)
634 {
635   Tao            tao = (Tao)ptr;
636   TAO_LCL        *lclP = (TAO_LCL*)tao->data;
637   PetscBool      set,pset,flag,pflag,symmetric;
638   PetscReal      cdotl;
639   PetscErrorCode ierr;
640 
641   PetscFunctionBegin;
642   ierr = TaoComputeObjectiveAndGradient(tao,X,f,G);CHKERRQ(ierr);
643   ierr = LCLScatter(lclP,G,lclP->GU,lclP->GV);CHKERRQ(ierr);
644   if (lclP->recompute_jacobian_flag) {
645     ierr = TaoComputeJacobianState(tao,X,tao->jacobian_state,tao->jacobian_state_pre,tao->jacobian_state_inv);CHKERRQ(ierr);
646     ierr = TaoComputeJacobianDesign(tao,X,tao->jacobian_design);CHKERRQ(ierr);
647   }
648   ierr = TaoComputeConstraints(tao,X, tao->constraints);CHKERRQ(ierr);
649   ierr = MatIsSymmetricKnown(tao->jacobian_state,&set,&flag);CHKERRQ(ierr);
650   if (tao->jacobian_state_pre) {
651     ierr = MatIsSymmetricKnown(tao->jacobian_state_pre,&pset,&pflag);CHKERRQ(ierr);
652   } else {
653     pset = pflag = PETSC_TRUE;
654   }
655   if (set && pset && flag && pflag) symmetric = PETSC_TRUE;
656   else symmetric = PETSC_FALSE;
657 
658   ierr = VecDot(lclP->lamda0, tao->constraints, &cdotl);CHKERRQ(ierr);
659   lclP->lgn = *f - cdotl;
660 
661   /* Gradient of Lagrangian GL = G - J' * lamda */
662   /*      WU = A' * WL
663           WV = B' * WL */
664   if (symmetric) {
665     ierr = MatMult(tao->jacobian_state,lclP->lamda0,lclP->GL_U);CHKERRQ(ierr);
666   } else {
667     ierr = MatMultTranspose(tao->jacobian_state,lclP->lamda0,lclP->GL_U);CHKERRQ(ierr);
668   }
669   ierr = MatMultTranspose(tao->jacobian_design,lclP->lamda0,lclP->GL_V);CHKERRQ(ierr);
670   ierr = VecScale(lclP->GL_U,-1.0);CHKERRQ(ierr);
671   ierr = VecScale(lclP->GL_V,-1.0);CHKERRQ(ierr);
672   ierr = VecAXPY(lclP->GL_U,1.0,lclP->GU);CHKERRQ(ierr);
673   ierr = VecAXPY(lclP->GL_V,1.0,lclP->GV);CHKERRQ(ierr);
674   ierr = LCLGather(lclP,lclP->GL_U,lclP->GL_V,lclP->GL);CHKERRQ(ierr);
675 
676   f[0] = lclP->lgn;
677   ierr = VecCopy(lclP->GL,G);CHKERRQ(ierr);
678   PetscFunctionReturn(0);
679 }
680 
681 static PetscErrorCode LCLComputeAugmentedLagrangianAndGradient(TaoLineSearch ls, Vec X, PetscReal *f, Vec G, void *ptr)
682 {
683   Tao            tao = (Tao)ptr;
684   TAO_LCL        *lclP = (TAO_LCL*)tao->data;
685   PetscReal      con2;
686   PetscBool      flag,pflag,set,pset,symmetric;
687   PetscErrorCode ierr;
688 
689   PetscFunctionBegin;
690   ierr = LCLComputeLagrangianAndGradient(tao->linesearch,X,f,G,tao);CHKERRQ(ierr);
691   ierr = LCLScatter(lclP,G,lclP->GL_U,lclP->GL_V);CHKERRQ(ierr);
692   ierr = VecDot(tao->constraints,tao->constraints,&con2);CHKERRQ(ierr);
693   lclP->aug = lclP->lgn + 0.5*lclP->rho*con2;
694 
695   /* Gradient of Aug. Lagrangian GAugL = GL + rho * J' c */
696   /*      WU = A' * c
697           WV = B' * c */
698   ierr = MatIsSymmetricKnown(tao->jacobian_state,&set,&flag);CHKERRQ(ierr);
699   if (tao->jacobian_state_pre) {
700     ierr = MatIsSymmetricKnown(tao->jacobian_state_pre,&pset,&pflag);CHKERRQ(ierr);
701   } else {
702     pset = pflag = PETSC_TRUE;
703   }
704   if (set && pset && flag && pflag) symmetric = PETSC_TRUE;
705   else symmetric = PETSC_FALSE;
706 
707   if (symmetric) {
708     ierr = MatMult(tao->jacobian_state,tao->constraints,lclP->GAugL_U);CHKERRQ(ierr);
709   } else {
710     ierr = MatMultTranspose(tao->jacobian_state,tao->constraints,lclP->GAugL_U);CHKERRQ(ierr);
711   }
712 
713   ierr = MatMultTranspose(tao->jacobian_design,tao->constraints,lclP->GAugL_V);CHKERRQ(ierr);
714   ierr = VecAYPX(lclP->GAugL_U,lclP->rho,lclP->GL_U);CHKERRQ(ierr);
715   ierr = VecAYPX(lclP->GAugL_V,lclP->rho,lclP->GL_V);CHKERRQ(ierr);
716   ierr = LCLGather(lclP,lclP->GAugL_U,lclP->GAugL_V,lclP->GAugL);CHKERRQ(ierr);
717 
718   f[0] = lclP->aug;
719   ierr = VecCopy(lclP->GAugL,G);CHKERRQ(ierr);
720   PetscFunctionReturn(0);
721 }
722 
723 PetscErrorCode LCLGather(TAO_LCL *lclP, Vec u, Vec v, Vec x)
724 {
725   PetscErrorCode ierr;
726   PetscFunctionBegin;
727   ierr = VecScatterBegin(lclP->state_scatter, u, x, INSERT_VALUES, SCATTER_REVERSE);CHKERRQ(ierr);
728   ierr = VecScatterEnd(lclP->state_scatter, u, x, INSERT_VALUES, SCATTER_REVERSE);CHKERRQ(ierr);
729   ierr = VecScatterBegin(lclP->design_scatter, v, x, INSERT_VALUES, SCATTER_REVERSE);CHKERRQ(ierr);
730   ierr = VecScatterEnd(lclP->design_scatter, v, x, INSERT_VALUES, SCATTER_REVERSE);CHKERRQ(ierr);
731   PetscFunctionReturn(0);
732 
733 }
734 PetscErrorCode LCLScatter(TAO_LCL *lclP, Vec x, Vec u, Vec v)
735 {
736   PetscErrorCode ierr;
737   PetscFunctionBegin;
738   ierr = VecScatterBegin(lclP->state_scatter, x, u, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr);
739   ierr = VecScatterEnd(lclP->state_scatter, x, u, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr);
740   ierr = VecScatterBegin(lclP->design_scatter, x, v, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr);
741   ierr = VecScatterEnd(lclP->design_scatter, x, v, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr);
742   PetscFunctionReturn(0);
743 
744 }
745