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