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