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