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