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