xref: /petsc/src/tao/pde_constrained/impls/lcl/lcl.c (revision 65f8aed5f7eaa1e2ef2ddeffe666264e0669c876)
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(PETSC_COMM_WORLD,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);
163     ierr = VecGetSize(lclP->V,&sizeV);
164     ierr = PetscPrintf(PETSC_COMM_WORLD,"size(U)=%D, size(V)=%D\n",sizeU,sizeV);
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     tao->ksp_its=0;
252     /* Compute a descent direction for the linearly constrained subproblem
253        minimize f(u+du, v+dv)
254        s.t. A(u0,v0)du + B(u0,v0)dv = -g(u0,v0) */
255 
256     /* Store the points around the linearization */
257     ierr = VecCopy(lclP->U, lclP->U0);CHKERRQ(ierr);
258     ierr = VecCopy(lclP->V, lclP->V0);CHKERRQ(ierr);
259     ierr = VecCopy(lclP->GU,lclP->GU0);CHKERRQ(ierr);
260     ierr = VecCopy(lclP->GV,lclP->GV0);CHKERRQ(ierr);
261     ierr = VecCopy(lclP->GAugL_U,lclP->GAugL_U0);CHKERRQ(ierr);
262     ierr = VecCopy(lclP->GAugL_V,lclP->GAugL_V0);CHKERRQ(ierr);
263     ierr = VecCopy(lclP->GL_U,lclP->GL_U0);CHKERRQ(ierr);
264     ierr = VecCopy(lclP->GL_V,lclP->GL_V0);CHKERRQ(ierr);
265 
266     lclP->aug0 = lclP->aug;
267     lclP->lgn0 = lclP->lgn;
268 
269     /* Given the design variables, we need to project the current iterate
270        onto the linearized constraint.  We choose to fix the design variables
271        and solve the linear system for the state variables.  The resulting
272        point is the Newton direction */
273 
274     /* Solve r = A\con */
275     lclP->solve_type = LCL_FORWARD1;
276     ierr = VecSet(lclP->r,0.0);CHKERRQ(ierr); /*  Initial guess in CG */
277 
278     if (tao->jacobian_state_inv) {
279       ierr = MatMult(tao->jacobian_state_inv, tao->constraints, lclP->r);CHKERRQ(ierr);
280     } else {
281       ierr = KSPSetOperators(tao->ksp, tao->jacobian_state, tao->jacobian_state_pre);CHKERRQ(ierr);
282       ierr = KSPSolve(tao->ksp, tao->constraints,  lclP->r);CHKERRQ(ierr);
283       ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr);
284       tao->ksp_its+=its;
285       tao->ksp_tot_its+=tao->ksp_its;
286     }
287 
288     /* Set design step direction dv to zero */
289     ierr = VecSet(lclP->s, 0.0);CHKERRQ(ierr);
290 
291     /*
292        Check sufficient descent for constraint merit function .5*||con||^2
293        con' Ak r >= eps1 ||r||^(2+eps2)
294     */
295 
296     /* Compute WU= Ak' * con */
297     if (symmetric)  {
298       ierr = MatMult(tao->jacobian_state,tao->constraints,lclP->WU);CHKERRQ(ierr);
299     } else {
300       ierr = MatMultTranspose(tao->jacobian_state,tao->constraints,lclP->WU);CHKERRQ(ierr);
301     }
302     /* Compute r * Ak' * con */
303     ierr = VecDot(lclP->r,lclP->WU,&rWU);CHKERRQ(ierr);
304 
305     /* compute ||r||^(2+eps2) */
306     ierr = VecNorm(lclP->r,NORM_2,&r2);CHKERRQ(ierr);
307     r2 = PetscPowScalar(r2,2.0+lclP->eps2);
308     adec = lclP->eps1 * r2;
309 
310     if (rWU < adec) {
311       ierr = PetscInfo(tao,"Newton direction not descent for constraint, feasibility phase required\n");CHKERRQ(ierr);
312       if (lclP->verbose) {
313         ierr = PetscPrintf(PETSC_COMM_WORLD,"Newton direction not descent for constraint: %g -- using steepest descent\n",(double)descent);CHKERRQ(ierr);
314       }
315 
316       ierr = PetscInfo(tao,"Using steepest descent direction instead.\n");CHKERRQ(ierr);
317       ierr = VecSet(lclP->r,0.0);CHKERRQ(ierr);
318       ierr = VecAXPY(lclP->r,-1.0,lclP->WU);CHKERRQ(ierr);
319       ierr = VecDot(lclP->r,lclP->r,&rWU);CHKERRQ(ierr);
320       ierr = VecNorm(lclP->r,NORM_2,&r2);CHKERRQ(ierr);
321       r2 = PetscPowScalar(r2,2.0+lclP->eps2);
322       ierr = VecDot(lclP->r,lclP->GAugL_U,&descent);CHKERRQ(ierr);
323       adec = lclP->eps1 * r2;
324     }
325 
326 
327     /*
328        Check descent for aug. lagrangian
329        r' (GUk - Ak'*yk - rho*Ak'*con) <= -eps1 ||r||^(2+eps2)
330           GL_U = GUk - Ak'*yk
331           WU   = Ak'*con
332                                          adec=eps1||r||^(2+eps2)
333 
334        ==>
335        Check r'GL_U - rho*r'WU <= adec
336     */
337 
338     ierr = VecDot(lclP->r,lclP->GL_U,&rGL_U);CHKERRQ(ierr);
339     aldescent =  rGL_U - lclP->rho*rWU;
340     if (aldescent > -adec) {
341       if (lclP->verbose) {
342         ierr = PetscPrintf(PETSC_COMM_WORLD," Newton direction not descent for augmented Lagrangian: %g",(double)aldescent);CHKERRQ(ierr);
343       }
344       ierr = PetscInfo1(tao,"Newton direction not descent for augmented Lagrangian: %g\n",(double)aldescent);CHKERRQ(ierr);
345       lclP->rho =  (rGL_U - adec)/rWU;
346       if (lclP->rho > lclP->rhomax) {
347         lclP->rho = lclP->rhomax;
348         SETERRQ1(PETSC_COMM_WORLD,0,"rho=%g > rhomax, case not implemented.  Increase rhomax (-tao_lcl_rhomax)",(double)lclP->rho);
349       }
350       if (lclP->verbose) {
351         ierr = PetscPrintf(PETSC_COMM_WORLD,"  Increasing penalty parameter to %g\n",(double)lclP->rho);CHKERRQ(ierr);
352       }
353       ierr = PetscInfo1(tao,"  Increasing penalty parameter to %g\n",(double)lclP->rho);CHKERRQ(ierr);
354     }
355 
356     ierr = LCLComputeAugmentedLagrangianAndGradient(tao->linesearch,tao->solution,&lclP->aug,lclP->GAugL,tao);CHKERRQ(ierr);
357     ierr = LCLScatter(lclP,lclP->GL,lclP->GL_U,lclP->GL_V);CHKERRQ(ierr);
358     ierr = LCLScatter(lclP,lclP->GAugL,lclP->GAugL_U,lclP->GAugL_V);CHKERRQ(ierr);
359 
360     /* We now minimize the augmented Lagrangian along the Newton direction */
361     ierr = VecScale(lclP->r,-1.0);CHKERRQ(ierr);
362     ierr = LCLGather(lclP, lclP->r,lclP->s,tao->stepdirection);CHKERRQ(ierr);
363     ierr = VecScale(lclP->r,-1.0);CHKERRQ(ierr);
364     ierr = LCLGather(lclP, lclP->GAugL_U0, lclP->GAugL_V0, lclP->GAugL);CHKERRQ(ierr);
365     ierr = LCLGather(lclP, lclP->U0,lclP->V0,lclP->X0);CHKERRQ(ierr);
366 
367     lclP->recompute_jacobian_flag = PETSC_TRUE;
368 
369     ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);CHKERRQ(ierr);
370     ierr = TaoLineSearchSetType(tao->linesearch, TAOLINESEARCHMT);CHKERRQ(ierr);
371     ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
372     ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &lclP->aug, lclP->GAugL, tao->stepdirection, &step, &ls_reason);CHKERRQ(ierr);
373     if (lclP->verbose) {
374       ierr = PetscPrintf(PETSC_COMM_WORLD,"Steplength = %g\n",(double)step);CHKERRQ(ierr);
375     }
376 
377     ierr = LCLScatter(lclP,tao->solution,lclP->U,lclP->V);CHKERRQ(ierr);
378     ierr = TaoComputeObjectiveAndGradient(tao,tao->solution,&f,tao->gradient);CHKERRQ(ierr);
379     ierr = LCLScatter(lclP,tao->gradient,lclP->GU,lclP->GV);CHKERRQ(ierr);
380 
381     ierr = LCLScatter(lclP,lclP->GAugL,lclP->GAugL_U,lclP->GAugL_V);CHKERRQ(ierr);
382 
383     /* Check convergence */
384     ierr = VecNorm(lclP->GAugL, NORM_2, &mnorm);CHKERRQ(ierr);
385     ierr = VecNorm(tao->constraints, NORM_2, &cnorm);CHKERRQ(ierr);
386     ierr = TaoLogConvergenceHistory(tao,f,mnorm,cnorm,tao->ksp_its);CHKERRQ(ierr);
387     ierr = TaoMonitor(tao,tao->niter,f,mnorm,cnorm,step);CHKERRQ(ierr);
388     ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
389     if (tao->reason != TAO_CONTINUE_ITERATING){
390       break;
391     }
392 
393     /* TODO: use a heuristic to choose how many iterations should be performed within phase 2 */
394     for (phase2_iter=0; phase2_iter<lclP->phase2_niter; phase2_iter++){
395       /* We now minimize the objective function starting from the fraction of
396          the Newton point accepted by applying one step of a reduced-space
397          method.  The optimization problem is:
398 
399          minimize f(u+du, v+dv)
400          s. t.    A(u0,v0)du + B(u0,v0)du = -alpha g(u0,v0)
401 
402          In particular, we have that
403          du = -inv(A)*(Bdv + alpha g) */
404 
405       ierr = TaoComputeJacobianState(tao,lclP->X0,tao->jacobian_state,tao->jacobian_state_pre,tao->jacobian_state_inv);CHKERRQ(ierr);
406       ierr = TaoComputeJacobianDesign(tao,lclP->X0,tao->jacobian_design);CHKERRQ(ierr);
407 
408       /* Store V and constraints */
409       ierr = VecCopy(lclP->V, lclP->V1);CHKERRQ(ierr);
410       ierr = VecCopy(tao->constraints,lclP->con1);CHKERRQ(ierr);
411 
412       /* Compute multipliers */
413       /* p1 */
414       ierr = VecSet(lclP->lamda,0.0);CHKERRQ(ierr); /*  Initial guess in CG */
415       lclP->solve_type = LCL_ADJOINT1;
416       ierr = MatIsSymmetricKnown(tao->jacobian_state,&set,&flag);CHKERRQ(ierr);
417       if (tao->jacobian_state_pre) {
418         ierr = MatIsSymmetricKnown(tao->jacobian_state_pre,&pset,&pflag);CHKERRQ(ierr);
419       } else {
420         pset = pflag = PETSC_TRUE;
421       }
422       if (set && pset && flag && pflag) symmetric = PETSC_TRUE;
423       else symmetric = PETSC_FALSE;
424 
425       if (tao->jacobian_state_inv) {
426         if (symmetric) {
427           ierr = MatMult(tao->jacobian_state_inv, lclP->GAugL_U, lclP->lamda);CHKERRQ(ierr);
428         } else {
429           ierr = MatMultTranspose(tao->jacobian_state_inv, lclP->GAugL_U, lclP->lamda);CHKERRQ(ierr);
430         }
431       } else {
432         if (symmetric) {
433           ierr = KSPSolve(tao->ksp, lclP->GAugL_U,  lclP->lamda);CHKERRQ(ierr);
434         } else {
435           ierr = KSPSolveTranspose(tao->ksp, lclP->GAugL_U,  lclP->lamda);CHKERRQ(ierr);
436         }
437         ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr);
438         tao->ksp_its+=its;
439         tao->ksp_tot_its+=its;
440       }
441       ierr = MatMultTranspose(tao->jacobian_design,lclP->lamda,lclP->g1);CHKERRQ(ierr);
442       ierr = VecAXPY(lclP->g1,-1.0,lclP->GAugL_V);CHKERRQ(ierr);
443 
444       /* Compute the limited-memory quasi-newton direction */
445       if (tao->niter > 0) {
446         ierr = MatSolve(lclP->R,lclP->g1,lclP->s);CHKERRQ(ierr);
447         ierr = VecDot(lclP->s,lclP->g1,&descent);CHKERRQ(ierr);
448         if (descent <= 0) {
449           if (lclP->verbose) {
450             ierr = PetscPrintf(PETSC_COMM_WORLD,"Reduced-space direction not descent: %g\n",(double)descent);CHKERRQ(ierr);
451           }
452           ierr = VecCopy(lclP->g1,lclP->s);CHKERRQ(ierr);
453         }
454       } else {
455         ierr = VecCopy(lclP->g1,lclP->s);CHKERRQ(ierr);
456       }
457       ierr = VecScale(lclP->g1,-1.0);CHKERRQ(ierr);
458 
459       /* Recover the full space direction */
460       ierr = MatMult(tao->jacobian_design,lclP->s,lclP->WU);CHKERRQ(ierr);
461       /* ierr = VecSet(lclP->r,0.0);CHKERRQ(ierr); */ /*  Initial guess in CG */
462       lclP->solve_type = LCL_FORWARD2;
463       if (tao->jacobian_state_inv) {
464         ierr = MatMult(tao->jacobian_state_inv,lclP->WU,lclP->r);CHKERRQ(ierr);
465       } else {
466         ierr = KSPSolve(tao->ksp, lclP->WU, lclP->r);CHKERRQ(ierr);
467         ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr);
468         tao->ksp_its += its;
469         tao->ksp_tot_its+=its;
470       }
471 
472       /* We now minimize the augmented Lagrangian along the direction -r,s */
473       ierr = VecScale(lclP->r, -1.0);CHKERRQ(ierr);
474       ierr = LCLGather(lclP,lclP->r,lclP->s,tao->stepdirection);CHKERRQ(ierr);
475       ierr = VecScale(lclP->r, -1.0);CHKERRQ(ierr);
476       lclP->recompute_jacobian_flag = PETSC_TRUE;
477 
478       ierr = TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);CHKERRQ(ierr);
479       ierr = TaoLineSearchSetType(tao->linesearch,TAOLINESEARCHMT);CHKERRQ(ierr);
480       ierr = TaoLineSearchSetFromOptions(tao->linesearch);CHKERRQ(ierr);
481       ierr = TaoLineSearchApply(tao->linesearch, tao->solution, &lclP->aug, lclP->GAugL, tao->stepdirection,&step,&ls_reason);CHKERRQ(ierr);
482       if (lclP->verbose){
483         ierr = PetscPrintf(PETSC_COMM_WORLD,"Reduced-space steplength =  %g\n",(double)step);CHKERRQ(ierr);
484       }
485 
486       ierr = LCLScatter(lclP,tao->solution,lclP->U,lclP->V);CHKERRQ(ierr);
487       ierr = LCLScatter(lclP,lclP->GL,lclP->GL_U,lclP->GL_V);CHKERRQ(ierr);
488       ierr = LCLScatter(lclP,lclP->GAugL,lclP->GAugL_U,lclP->GAugL_V);CHKERRQ(ierr);
489       ierr = TaoComputeObjectiveAndGradient(tao,tao->solution,&f,tao->gradient);CHKERRQ(ierr);
490       ierr = LCLScatter(lclP,tao->gradient,lclP->GU,lclP->GV);CHKERRQ(ierr);
491 
492       /* Compute the reduced gradient at the new point */
493 
494       ierr = TaoComputeJacobianState(tao,lclP->X0,tao->jacobian_state,tao->jacobian_state_pre,tao->jacobian_state_inv);CHKERRQ(ierr);
495       ierr = TaoComputeJacobianDesign(tao,lclP->X0,tao->jacobian_design);CHKERRQ(ierr);
496 
497       /* p2 */
498       /* Compute multipliers, using lamda-rho*con as an initial guess in PCG */
499       if (phase2_iter==0){
500         ierr = VecWAXPY(lclP->lamda,-lclP->rho,lclP->con1,lclP->lamda0);CHKERRQ(ierr);
501       } else {
502         ierr = VecAXPY(lclP->lamda,-lclP->rho,tao->constraints);CHKERRQ(ierr);
503       }
504 
505       ierr = MatIsSymmetricKnown(tao->jacobian_state,&set,&flag);CHKERRQ(ierr);
506       if (tao->jacobian_state_pre) {
507         ierr = MatIsSymmetricKnown(tao->jacobian_state_pre,&pset,&pflag);CHKERRQ(ierr);
508       } else {
509         pset = pflag = PETSC_TRUE;
510       }
511       if (set && pset && flag && pflag) symmetric = PETSC_TRUE;
512       else symmetric = PETSC_FALSE;
513 
514       lclP->solve_type = LCL_ADJOINT2;
515       if (tao->jacobian_state_inv) {
516         if (symmetric) {
517           ierr = MatMult(tao->jacobian_state_inv, lclP->GU, lclP->lamda);CHKERRQ(ierr);
518         } else {
519           ierr = MatMultTranspose(tao->jacobian_state_inv, lclP->GU, lclP->lamda);CHKERRQ(ierr);
520         }
521       } else {
522         if (symmetric) {
523           ierr = KSPSolve(tao->ksp, lclP->GU,  lclP->lamda);CHKERRQ(ierr);
524         } else {
525           ierr = KSPSolveTranspose(tao->ksp, lclP->GU,  lclP->lamda);CHKERRQ(ierr);
526         }
527         ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr);
528         tao->ksp_its += its;
529         tao->ksp_tot_its += its;
530       }
531 
532       ierr = MatMultTranspose(tao->jacobian_design,lclP->lamda,lclP->g2);CHKERRQ(ierr);
533       ierr = VecAXPY(lclP->g2,-1.0,lclP->GV);CHKERRQ(ierr);
534 
535       ierr = VecScale(lclP->g2,-1.0);CHKERRQ(ierr);
536 
537       /* Update the quasi-newton approximation */
538       ierr = MatLMVMUpdate(lclP->R,lclP->V,lclP->g2);CHKERRQ(ierr);
539       /* 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 */
540 
541     }
542 
543     ierr = VecCopy(lclP->lamda,lclP->lamda0);CHKERRQ(ierr);
544 
545     /* Evaluate Function, Gradient, Constraints, and Jacobian */
546     ierr = TaoComputeObjectiveAndGradient(tao,tao->solution,&f,tao->gradient);CHKERRQ(ierr);
547     ierr = LCLScatter(lclP,tao->solution,lclP->U,lclP->V);CHKERRQ(ierr);
548     ierr = LCLScatter(lclP,tao->gradient,lclP->GU,lclP->GV);CHKERRQ(ierr);
549 
550     ierr = TaoComputeJacobianState(tao,tao->solution,tao->jacobian_state,tao->jacobian_state_pre,tao->jacobian_state_inv);CHKERRQ(ierr);
551     ierr = TaoComputeJacobianDesign(tao,tao->solution,tao->jacobian_design);CHKERRQ(ierr);
552     ierr = TaoComputeConstraints(tao,tao->solution, tao->constraints);CHKERRQ(ierr);
553 
554     ierr = LCLComputeAugmentedLagrangianAndGradient(tao->linesearch,tao->solution,&lclP->aug,lclP->GAugL,tao);CHKERRQ(ierr);
555 
556     ierr = VecNorm(lclP->GAugL, NORM_2, &mnorm);CHKERRQ(ierr);
557 
558     /* Evaluate constraint norm */
559     ierr = VecNorm(tao->constraints, NORM_2, &cnorm);CHKERRQ(ierr);
560 
561     /* Monitor convergence */
562     tao->niter++;
563     ierr = TaoLogConvergenceHistory(tao,f,mnorm,cnorm,tao->ksp_its);CHKERRQ(ierr);
564     ierr = TaoMonitor(tao,tao->niter,f,mnorm,cnorm,step);CHKERRQ(ierr);
565     ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
566   }
567   PetscFunctionReturn(0);
568 }
569 
570 /*MC
571  TAOLCL - linearly constrained lagrangian method for pde-constrained optimization
572 
573 + -tao_lcl_eps1 - epsilon 1 tolerance
574 . -tao_lcl_eps2","epsilon 2 tolerance","",lclP->eps2,&lclP->eps2,NULL);CHKERRQ(ierr);
575 . -tao_lcl_rho0","init value for rho","",lclP->rho0,&lclP->rho0,NULL);CHKERRQ(ierr);
576 . -tao_lcl_rhomax","max value for rho","",lclP->rhomax,&lclP->rhomax,NULL);CHKERRQ(ierr);
577 . -tao_lcl_phase2_niter - Number of phase 2 iterations in LCL algorithm
578 . -tao_lcl_verbose - Print verbose output if True
579 . -tao_lcl_tola - Tolerance for first forward solve
580 . -tao_lcl_tolb - Tolerance for first adjoint solve
581 . -tao_lcl_tolc - Tolerance for second forward solve
582 - -tao_lcl_told - Tolerance for second adjoint solve
583 
584   Level: beginner
585 M*/
586 PETSC_EXTERN PetscErrorCode TaoCreate_LCL(Tao tao)
587 {
588   TAO_LCL        *lclP;
589   PetscErrorCode ierr;
590   const char     *morethuente_type = TAOLINESEARCHMT;
591 
592   PetscFunctionBegin;
593   tao->ops->setup = TaoSetup_LCL;
594   tao->ops->solve = TaoSolve_LCL;
595   tao->ops->view = TaoView_LCL;
596   tao->ops->setfromoptions = TaoSetFromOptions_LCL;
597   tao->ops->destroy = TaoDestroy_LCL;
598   ierr = PetscNewLog(tao,&lclP);CHKERRQ(ierr);
599   tao->data = (void*)lclP;
600 
601   /* Override default settings (unless already changed) */
602   if (!tao->max_it_changed) tao->max_it = 200;
603   if (!tao->catol_changed) tao->catol = 1.0e-4;
604   if (!tao->gatol_changed) tao->gttol = 1.0e-4;
605   if (!tao->grtol_changed) tao->gttol = 1.0e-4;
606   if (!tao->gttol_changed) tao->gttol = 1.0e-4;
607   lclP->rho0 = 1.0e-4;
608   lclP->rhomax=1e5;
609   lclP->eps1 = 1.0e-8;
610   lclP->eps2 = 0.0;
611   lclP->solve_type=2;
612   lclP->tau[0] = lclP->tau[1] = lclP->tau[2] = lclP->tau[3] = 1.0e-4;
613   ierr = TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);CHKERRQ(ierr);
614   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);CHKERRQ(ierr);
615   ierr = TaoLineSearchSetType(tao->linesearch, morethuente_type);CHKERRQ(ierr);
616   ierr = TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);CHKERRQ(ierr);
617 
618   ierr = TaoLineSearchSetObjectiveAndGradientRoutine(tao->linesearch,LCLComputeAugmentedLagrangianAndGradient, tao);CHKERRQ(ierr);
619   ierr = KSPCreate(((PetscObject)tao)->comm,&tao->ksp);CHKERRQ(ierr);
620   ierr = PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1);CHKERRQ(ierr);
621   ierr = KSPSetOptionsPrefix(tao->ksp, tao->hdr.prefix);CHKERRQ(ierr);
622   ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
623 
624   ierr = MatCreate(((PetscObject)tao)->comm, &lclP->R);CHKERRQ(ierr);
625   ierr = MatSetType(lclP->R, MATLMVMBFGS);CHKERRQ(ierr);
626   PetscFunctionReturn(0);
627 }
628 
629 static PetscErrorCode LCLComputeLagrangianAndGradient(TaoLineSearch ls, Vec X, PetscReal *f, Vec G, void *ptr)
630 {
631   Tao            tao = (Tao)ptr;
632   TAO_LCL        *lclP = (TAO_LCL*)tao->data;
633   PetscBool      set,pset,flag,pflag,symmetric;
634   PetscReal      cdotl;
635   PetscErrorCode ierr;
636 
637   PetscFunctionBegin;
638   ierr = TaoComputeObjectiveAndGradient(tao,X,f,G);CHKERRQ(ierr);
639   ierr = LCLScatter(lclP,G,lclP->GU,lclP->GV);CHKERRQ(ierr);
640   if (lclP->recompute_jacobian_flag) {
641     ierr = TaoComputeJacobianState(tao,X,tao->jacobian_state,tao->jacobian_state_pre,tao->jacobian_state_inv);CHKERRQ(ierr);
642     ierr = TaoComputeJacobianDesign(tao,X,tao->jacobian_design);CHKERRQ(ierr);
643   }
644   ierr = TaoComputeConstraints(tao,X, tao->constraints);CHKERRQ(ierr);
645   ierr = MatIsSymmetricKnown(tao->jacobian_state,&set,&flag);CHKERRQ(ierr);
646   if (tao->jacobian_state_pre) {
647     ierr = MatIsSymmetricKnown(tao->jacobian_state_pre,&pset,&pflag);CHKERRQ(ierr);
648   } else {
649     pset = pflag = PETSC_TRUE;
650   }
651   if (set && pset && flag && pflag) symmetric = PETSC_TRUE;
652   else symmetric = PETSC_FALSE;
653 
654   ierr = VecDot(lclP->lamda0, tao->constraints, &cdotl);CHKERRQ(ierr);
655   lclP->lgn = *f - cdotl;
656 
657   /* Gradient of Lagrangian GL = G - J' * lamda */
658   /*      WU = A' * WL
659           WV = B' * WL */
660   if (symmetric) {
661     ierr = MatMult(tao->jacobian_state,lclP->lamda0,lclP->GL_U);CHKERRQ(ierr);
662   } else {
663     ierr = MatMultTranspose(tao->jacobian_state,lclP->lamda0,lclP->GL_U);CHKERRQ(ierr);
664   }
665   ierr = MatMultTranspose(tao->jacobian_design,lclP->lamda0,lclP->GL_V);CHKERRQ(ierr);
666   ierr = VecScale(lclP->GL_U,-1.0);CHKERRQ(ierr);
667   ierr = VecScale(lclP->GL_V,-1.0);CHKERRQ(ierr);
668   ierr = VecAXPY(lclP->GL_U,1.0,lclP->GU);CHKERRQ(ierr);
669   ierr = VecAXPY(lclP->GL_V,1.0,lclP->GV);CHKERRQ(ierr);
670   ierr = LCLGather(lclP,lclP->GL_U,lclP->GL_V,lclP->GL);CHKERRQ(ierr);
671 
672   f[0] = lclP->lgn;
673   ierr = VecCopy(lclP->GL,G);CHKERRQ(ierr);
674   PetscFunctionReturn(0);
675 }
676 
677 static PetscErrorCode LCLComputeAugmentedLagrangianAndGradient(TaoLineSearch ls, Vec X, PetscReal *f, Vec G, void *ptr)
678 {
679   Tao            tao = (Tao)ptr;
680   TAO_LCL        *lclP = (TAO_LCL*)tao->data;
681   PetscReal      con2;
682   PetscBool      flag,pflag,set,pset,symmetric;
683   PetscErrorCode ierr;
684 
685   PetscFunctionBegin;
686   ierr = LCLComputeLagrangianAndGradient(tao->linesearch,X,f,G,tao);CHKERRQ(ierr);
687   ierr = LCLScatter(lclP,G,lclP->GL_U,lclP->GL_V);CHKERRQ(ierr);
688   ierr = VecDot(tao->constraints,tao->constraints,&con2);CHKERRQ(ierr);
689   lclP->aug = lclP->lgn + 0.5*lclP->rho*con2;
690 
691   /* Gradient of Aug. Lagrangian GAugL = GL + rho * J' c */
692   /*      WU = A' * c
693           WV = B' * c */
694   ierr = MatIsSymmetricKnown(tao->jacobian_state,&set,&flag);CHKERRQ(ierr);
695   if (tao->jacobian_state_pre) {
696     ierr = MatIsSymmetricKnown(tao->jacobian_state_pre,&pset,&pflag);CHKERRQ(ierr);
697   } else {
698     pset = pflag = PETSC_TRUE;
699   }
700   if (set && pset && flag && pflag) symmetric = PETSC_TRUE;
701   else symmetric = PETSC_FALSE;
702 
703   if (symmetric) {
704     ierr = MatMult(tao->jacobian_state,tao->constraints,lclP->GAugL_U);CHKERRQ(ierr);
705   } else {
706     ierr = MatMultTranspose(tao->jacobian_state,tao->constraints,lclP->GAugL_U);CHKERRQ(ierr);
707   }
708 
709   ierr = MatMultTranspose(tao->jacobian_design,tao->constraints,lclP->GAugL_V);CHKERRQ(ierr);
710   ierr = VecAYPX(lclP->GAugL_U,lclP->rho,lclP->GL_U);CHKERRQ(ierr);
711   ierr = VecAYPX(lclP->GAugL_V,lclP->rho,lclP->GL_V);CHKERRQ(ierr);
712   ierr = LCLGather(lclP,lclP->GAugL_U,lclP->GAugL_V,lclP->GAugL);CHKERRQ(ierr);
713 
714   f[0] = lclP->aug;
715   ierr = VecCopy(lclP->GAugL,G);CHKERRQ(ierr);
716   PetscFunctionReturn(0);
717 }
718 
719 PetscErrorCode LCLGather(TAO_LCL *lclP, Vec u, Vec v, Vec x)
720 {
721   PetscErrorCode ierr;
722   PetscFunctionBegin;
723   ierr = VecScatterBegin(lclP->state_scatter, u, x, INSERT_VALUES, SCATTER_REVERSE);CHKERRQ(ierr);
724   ierr = VecScatterEnd(lclP->state_scatter, u, x, INSERT_VALUES, SCATTER_REVERSE);CHKERRQ(ierr);
725   ierr = VecScatterBegin(lclP->design_scatter, v, x, INSERT_VALUES, SCATTER_REVERSE);CHKERRQ(ierr);
726   ierr = VecScatterEnd(lclP->design_scatter, v, x, INSERT_VALUES, SCATTER_REVERSE);CHKERRQ(ierr);
727   PetscFunctionReturn(0);
728 
729 }
730 PetscErrorCode LCLScatter(TAO_LCL *lclP, Vec x, Vec u, Vec v)
731 {
732   PetscErrorCode ierr;
733   PetscFunctionBegin;
734   ierr = VecScatterBegin(lclP->state_scatter, x, u, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr);
735   ierr = VecScatterEnd(lclP->state_scatter, x, u, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr);
736   ierr = VecScatterBegin(lclP->design_scatter, x, v, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr);
737   ierr = VecScatterEnd(lclP->design_scatter, x, v, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr);
738   PetscFunctionReturn(0);
739 
740 }
741