Actual source code: lcl.c

petsc-3.10.5 2019-03-28
Report Typos and Errors
  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);

  7: static PetscErrorCode TaoDestroy_LCL(Tao tao)
  8: {
  9:   TAO_LCL        *lclP = (TAO_LCL*)tao->data;

 13:   if (tao->setupcalled) {
 14:     MatDestroy(&lclP->R);
 15:     VecDestroy(&lclP->lamda);
 16:     VecDestroy(&lclP->lamda0);
 17:     VecDestroy(&lclP->WL);
 18:     VecDestroy(&lclP->W);
 19:     VecDestroy(&lclP->X0);
 20:     VecDestroy(&lclP->G0);
 21:     VecDestroy(&lclP->GL);
 22:     VecDestroy(&lclP->GAugL);
 23:     VecDestroy(&lclP->dbar);
 24:     VecDestroy(&lclP->U);
 25:     VecDestroy(&lclP->U0);
 26:     VecDestroy(&lclP->V);
 27:     VecDestroy(&lclP->V0);
 28:     VecDestroy(&lclP->V1);
 29:     VecDestroy(&lclP->GU);
 30:     VecDestroy(&lclP->GV);
 31:     VecDestroy(&lclP->GU0);
 32:     VecDestroy(&lclP->GV0);
 33:     VecDestroy(&lclP->GL_U);
 34:     VecDestroy(&lclP->GL_V);
 35:     VecDestroy(&lclP->GAugL_U);
 36:     VecDestroy(&lclP->GAugL_V);
 37:     VecDestroy(&lclP->GL_U0);
 38:     VecDestroy(&lclP->GL_V0);
 39:     VecDestroy(&lclP->GAugL_U0);
 40:     VecDestroy(&lclP->GAugL_V0);
 41:     VecDestroy(&lclP->DU);
 42:     VecDestroy(&lclP->DV);
 43:     VecDestroy(&lclP->WU);
 44:     VecDestroy(&lclP->WV);
 45:     VecDestroy(&lclP->g1);
 46:     VecDestroy(&lclP->g2);
 47:     VecDestroy(&lclP->con1);

 49:     VecDestroy(&lclP->r);
 50:     VecDestroy(&lclP->s);

 52:     ISDestroy(&tao->state_is);
 53:     ISDestroy(&tao->design_is);

 55:     VecScatterDestroy(&lclP->state_scatter);
 56:     VecScatterDestroy(&lclP->design_scatter);
 57:   }
 58:   MatDestroy(&lclP->R);
 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:   MatSetFromOptions(lclP->R);
 86:   return(0);
 87: }

 89: static PetscErrorCode TaoView_LCL(Tao tao, PetscViewer viewer)
 90: {
 91:   return 0;
 92: }

 94: static PetscErrorCode TaoSetup_LCL(Tao tao)
 95: {
 96:   TAO_LCL        *lclP = (TAO_LCL*)tao->data;
 97:   PetscInt       lo, hi, nlocalstate, nlocaldesign;
 99:   IS             is_state, is_design;

102:   if (!tao->state_is) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONGSTATE,"LCL Solver requires an initial state index set -- use TaoSetStateIS()");
103:   VecDuplicate(tao->solution, &tao->gradient);
104:   VecDuplicate(tao->solution, &tao->stepdirection);
105:   VecDuplicate(tao->solution, &lclP->W);
106:   VecDuplicate(tao->solution, &lclP->X0);
107:   VecDuplicate(tao->solution, &lclP->G0);
108:   VecDuplicate(tao->solution, &lclP->GL);
109:   VecDuplicate(tao->solution, &lclP->GAugL);

111:   VecDuplicate(tao->constraints, &lclP->lamda);
112:   VecDuplicate(tao->constraints, &lclP->WL);
113:   VecDuplicate(tao->constraints, &lclP->lamda0);
114:   VecDuplicate(tao->constraints, &lclP->con1);

116:   VecSet(lclP->lamda,0.0);

118:   VecGetSize(tao->solution, &lclP->n);
119:   VecGetSize(tao->constraints, &lclP->m);

121:   VecCreate(((PetscObject)tao)->comm,&lclP->U);
122:   VecCreate(((PetscObject)tao)->comm,&lclP->V);
123:   ISGetLocalSize(tao->state_is,&nlocalstate);
124:   ISGetLocalSize(tao->design_is,&nlocaldesign);
125:   VecSetSizes(lclP->U,nlocalstate,lclP->m);
126:   VecSetSizes(lclP->V,nlocaldesign,lclP->n-lclP->m);
127:   VecSetType(lclP->U,((PetscObject)(tao->solution))->type_name);
128:   VecSetType(lclP->V,((PetscObject)(tao->solution))->type_name);
129:   VecSetFromOptions(lclP->U);
130:   VecSetFromOptions(lclP->V);
131:   VecDuplicate(lclP->U,&lclP->DU);
132:   VecDuplicate(lclP->U,&lclP->U0);
133:   VecDuplicate(lclP->U,&lclP->GU);
134:   VecDuplicate(lclP->U,&lclP->GU0);
135:   VecDuplicate(lclP->U,&lclP->GAugL_U);
136:   VecDuplicate(lclP->U,&lclP->GL_U);
137:   VecDuplicate(lclP->U,&lclP->GAugL_U0);
138:   VecDuplicate(lclP->U,&lclP->GL_U0);
139:   VecDuplicate(lclP->U,&lclP->WU);
140:   VecDuplicate(lclP->U,&lclP->r);
141:   VecDuplicate(lclP->V,&lclP->V0);
142:   VecDuplicate(lclP->V,&lclP->V1);
143:   VecDuplicate(lclP->V,&lclP->DV);
144:   VecDuplicate(lclP->V,&lclP->s);
145:   VecDuplicate(lclP->V,&lclP->GV);
146:   VecDuplicate(lclP->V,&lclP->GV0);
147:   VecDuplicate(lclP->V,&lclP->dbar);
148:   VecDuplicate(lclP->V,&lclP->GAugL_V);
149:   VecDuplicate(lclP->V,&lclP->GL_V);
150:   VecDuplicate(lclP->V,&lclP->GAugL_V0);
151:   VecDuplicate(lclP->V,&lclP->GL_V0);
152:   VecDuplicate(lclP->V,&lclP->WV);
153:   VecDuplicate(lclP->V,&lclP->g1);
154:   VecDuplicate(lclP->V,&lclP->g2);

156:   /* create scatters for state, design subvecs */
157:   VecGetOwnershipRange(lclP->U,&lo,&hi);
158:   ISCreateStride(((PetscObject)lclP->U)->comm,hi-lo,lo,1,&is_state);
159:   VecGetOwnershipRange(lclP->V,&lo,&hi);
160:   if (0) {
161:     PetscInt sizeU,sizeV;
162:     VecGetSize(lclP->U,&sizeU);
163:     VecGetSize(lclP->V,&sizeV);
164:     PetscPrintf(PETSC_COMM_WORLD,"size(U)=%D, size(V)=%D\n",sizeU,sizeV);
165:   }
166:   ISCreateStride(((PetscObject)lclP->V)->comm,hi-lo,lo,1,&is_design);
167:   VecScatterCreate(tao->solution,tao->state_is,lclP->U,is_state,&lclP->state_scatter);
168:   VecScatterCreate(tao->solution,tao->design_is,lclP->V,is_design,&lclP->design_scatter);
169:   ISDestroy(&is_state);
170:   ISDestroy(&is_design);
171:   return(0);
172: }

174: static PetscErrorCode TaoSolve_LCL(Tao tao)
175: {
176:   TAO_LCL                      *lclP = (TAO_LCL*)tao->data;
177:   PetscInt                     phase2_iter,nlocal,its;
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:   MatSetSizes(lclP->R, nlocal, nlocal, lclP->n-lclP->m, lclP->n-lclP->m);
190:   MatLMVMAllocate(lclP->R,lclP->V,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:   tao->reason = TAO_CONTINUE_ITERATING;
246:   TaoLogConvergenceHistory(tao,f,mnorm,cnorm,tao->ksp_its);
247:   TaoMonitor(tao,tao->niter,f,mnorm,cnorm,step);
248:   (*tao->ops->convergencetest)(tao,tao->cnvP);

250:   while (tao->reason == TAO_CONTINUE_ITERATING) {
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) */

256:     /* Store the points around the linearization */
257:     VecCopy(lclP->U, lclP->U0);
258:     VecCopy(lclP->V, lclP->V0);
259:     VecCopy(lclP->GU,lclP->GU0);
260:     VecCopy(lclP->GV,lclP->GV0);
261:     VecCopy(lclP->GAugL_U,lclP->GAugL_U0);
262:     VecCopy(lclP->GAugL_V,lclP->GAugL_V0);
263:     VecCopy(lclP->GL_U,lclP->GL_U0);
264:     VecCopy(lclP->GL_V,lclP->GL_V0);

266:     lclP->aug0 = lclP->aug;
267:     lclP->lgn0 = lclP->lgn;

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 */

274:     /* Solve r = A\con */
275:     lclP->solve_type = LCL_FORWARD1;
276:     VecSet(lclP->r,0.0); /*  Initial guess in CG */

278:     if (tao->jacobian_state_inv) {
279:       MatMult(tao->jacobian_state_inv, tao->constraints, lclP->r);
280:     } else {
281:       KSPSetOperators(tao->ksp, tao->jacobian_state, tao->jacobian_state_pre);
282:       KSPSolve(tao->ksp, tao->constraints,  lclP->r);
283:       KSPGetIterationNumber(tao->ksp,&its);
284:       tao->ksp_its+=its;
285:       tao->ksp_tot_its+=tao->ksp_its;
286:     }

288:     /* Set design step direction dv to zero */
289:     VecSet(lclP->s, 0.0);

291:     /*
292:        Check sufficient descent for constraint merit function .5*||con||^2
293:        con' Ak r >= eps1 ||r||^(2+eps2)
294:     */

296:     /* Compute WU= Ak' * con */
297:     if (symmetric)  {
298:       MatMult(tao->jacobian_state,tao->constraints,lclP->WU);
299:     } else {
300:       MatMultTranspose(tao->jacobian_state,tao->constraints,lclP->WU);
301:     }
302:     /* Compute r * Ak' * con */
303:     VecDot(lclP->r,lclP->WU,&rWU);

305:     /* compute ||r||^(2+eps2) */
306:     VecNorm(lclP->r,NORM_2,&r2);
307:     r2 = PetscPowScalar(r2,2.0+lclP->eps2);
308:     adec = lclP->eps1 * r2;

310:     if (rWU < adec) {
311:       PetscInfo(tao,"Newton direction not descent for constraint, feasibility phase required\n");
312:       if (lclP->verbose) {
313:         PetscPrintf(PETSC_COMM_WORLD,"Newton direction not descent for constraint: %g -- using steepest descent\n",(double)descent);
314:       }

316:       PetscInfo(tao,"Using steepest descent direction instead.\n");
317:       VecSet(lclP->r,0.0);
318:       VecAXPY(lclP->r,-1.0,lclP->WU);
319:       VecDot(lclP->r,lclP->r,&rWU);
320:       VecNorm(lclP->r,NORM_2,&r2);
321:       r2 = PetscPowScalar(r2,2.0+lclP->eps2);
322:       VecDot(lclP->r,lclP->GAugL_U,&descent);
323:       adec = lclP->eps1 * r2;
324:     }


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)

334:        ==>
335:        Check r'GL_U - rho*r'WU <= adec
336:     */

338:     VecDot(lclP->r,lclP->GL_U,&rGL_U);
339:     aldescent =  rGL_U - lclP->rho*rWU;
340:     if (aldescent > -adec) {
341:       if (lclP->verbose) {
342:         PetscPrintf(PETSC_COMM_WORLD," Newton direction not descent for augmented Lagrangian: %g",(double)aldescent);
343:       }
344:       PetscInfo1(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:         SETERRQ1(PETSC_COMM_WORLD,0,"rho=%g > rhomax, case not implemented.  Increase rhomax (-tao_lcl_rhomax)",(double)lclP->rho);
349:       }
350:       if (lclP->verbose) {
351:         PetscPrintf(PETSC_COMM_WORLD,"  Increasing penalty parameter to %g\n",(double)lclP->rho);
352:       }
353:       PetscInfo1(tao,"  Increasing penalty parameter to %g\n",(double)lclP->rho);
354:     }

356:     LCLComputeAugmentedLagrangianAndGradient(tao->linesearch,tao->solution,&lclP->aug,lclP->GAugL,tao);
357:     LCLScatter(lclP,lclP->GL,lclP->GL_U,lclP->GL_V);
358:     LCLScatter(lclP,lclP->GAugL,lclP->GAugL_U,lclP->GAugL_V);

360:     /* We now minimize the augmented Lagrangian along the Newton direction */
361:     VecScale(lclP->r,-1.0);
362:     LCLGather(lclP, lclP->r,lclP->s,tao->stepdirection);
363:     VecScale(lclP->r,-1.0);
364:     LCLGather(lclP, lclP->GAugL_U0, lclP->GAugL_V0, lclP->GAugL);
365:     LCLGather(lclP, lclP->U0,lclP->V0,lclP->X0);

367:     lclP->recompute_jacobian_flag = PETSC_TRUE;

369:     TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);
370:     TaoLineSearchSetType(tao->linesearch, TAOLINESEARCHMT);
371:     TaoLineSearchSetFromOptions(tao->linesearch);
372:     TaoLineSearchApply(tao->linesearch, tao->solution, &lclP->aug, lclP->GAugL, tao->stepdirection, &step, &ls_reason);
373:     if (lclP->verbose) {
374:       PetscPrintf(PETSC_COMM_WORLD,"Steplength = %g\n",(double)step);
375:     }

377:     LCLScatter(lclP,tao->solution,lclP->U,lclP->V);
378:     TaoComputeObjectiveAndGradient(tao,tao->solution,&f,tao->gradient);
379:     LCLScatter(lclP,tao->gradient,lclP->GU,lclP->GV);

381:     LCLScatter(lclP,lclP->GAugL,lclP->GAugL_U,lclP->GAugL_V);

383:     /* Check convergence */
384:     VecNorm(lclP->GAugL, NORM_2, &mnorm);
385:     VecNorm(tao->constraints, NORM_2, &cnorm);
386:     TaoLogConvergenceHistory(tao,f,mnorm,cnorm,tao->ksp_its);
387:     TaoMonitor(tao,tao->niter,f,mnorm,cnorm,step);
388:     (*tao->ops->convergencetest)(tao,tao->cnvP);
389:     if (tao->reason != TAO_CONTINUE_ITERATING){
390:       break;
391:     }

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:

399:          minimize f(u+du, v+dv)
400:          s. t.    A(u0,v0)du + B(u0,v0)du = -alpha g(u0,v0)

402:          In particular, we have that
403:          du = -inv(A)*(Bdv + alpha g) */

405:       TaoComputeJacobianState(tao,lclP->X0,tao->jacobian_state,tao->jacobian_state_pre,tao->jacobian_state_inv);
406:       TaoComputeJacobianDesign(tao,lclP->X0,tao->jacobian_design);

408:       /* Store V and constraints */
409:       VecCopy(lclP->V, lclP->V1);
410:       VecCopy(tao->constraints,lclP->con1);

412:       /* Compute multipliers */
413:       /* p1 */
414:       VecSet(lclP->lamda,0.0); /*  Initial guess in CG */
415:       lclP->solve_type = LCL_ADJOINT1;
416:       MatIsSymmetricKnown(tao->jacobian_state,&set,&flag);
417:       if (tao->jacobian_state_pre) {
418:         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;

425:       if (tao->jacobian_state_inv) {
426:         if (symmetric) {
427:           MatMult(tao->jacobian_state_inv, lclP->GAugL_U, lclP->lamda);
428:         } else {
429:           MatMultTranspose(tao->jacobian_state_inv, lclP->GAugL_U, lclP->lamda);
430:         }
431:       } else {
432:         if (symmetric) {
433:           KSPSolve(tao->ksp, lclP->GAugL_U,  lclP->lamda);
434:         } else {
435:           KSPSolveTranspose(tao->ksp, lclP->GAugL_U,  lclP->lamda);
436:         }
437:         KSPGetIterationNumber(tao->ksp,&its);
438:         tao->ksp_its+=its;
439:         tao->ksp_tot_its+=its;
440:       }
441:       MatMultTranspose(tao->jacobian_design,lclP->lamda,lclP->g1);
442:       VecAXPY(lclP->g1,-1.0,lclP->GAugL_V);

444:       /* Compute the limited-memory quasi-newton direction */
445:       if (tao->niter > 0) {
446:         MatSolve(lclP->R,lclP->g1,lclP->s);
447:         VecDot(lclP->s,lclP->g1,&descent);
448:         if (descent <= 0) {
449:           if (lclP->verbose) {
450:             PetscPrintf(PETSC_COMM_WORLD,"Reduced-space direction not descent: %g\n",(double)descent);
451:           }
452:           VecCopy(lclP->g1,lclP->s);
453:         }
454:       } else {
455:         VecCopy(lclP->g1,lclP->s);
456:       }
457:       VecScale(lclP->g1,-1.0);

459:       /* Recover the full space direction */
460:       MatMult(tao->jacobian_design,lclP->s,lclP->WU);
461:       /* VecSet(lclP->r,0.0); */ /*  Initial guess in CG */
462:       lclP->solve_type = LCL_FORWARD2;
463:       if (tao->jacobian_state_inv) {
464:         MatMult(tao->jacobian_state_inv,lclP->WU,lclP->r);
465:       } else {
466:         KSPSolve(tao->ksp, lclP->WU, lclP->r);
467:         KSPGetIterationNumber(tao->ksp,&its);
468:         tao->ksp_its += its;
469:         tao->ksp_tot_its+=its;
470:       }

472:       /* We now minimize the augmented Lagrangian along the direction -r,s */
473:       VecScale(lclP->r, -1.0);
474:       LCLGather(lclP,lclP->r,lclP->s,tao->stepdirection);
475:       VecScale(lclP->r, -1.0);
476:       lclP->recompute_jacobian_flag = PETSC_TRUE;

478:       TaoLineSearchSetInitialStepLength(tao->linesearch,1.0);
479:       TaoLineSearchSetType(tao->linesearch,TAOLINESEARCHMT);
480:       TaoLineSearchSetFromOptions(tao->linesearch);
481:       TaoLineSearchApply(tao->linesearch, tao->solution, &lclP->aug, lclP->GAugL, tao->stepdirection,&step,&ls_reason);
482:       if (lclP->verbose){
483:         PetscPrintf(PETSC_COMM_WORLD,"Reduced-space steplength =  %g\n",(double)step);
484:       }

486:       LCLScatter(lclP,tao->solution,lclP->U,lclP->V);
487:       LCLScatter(lclP,lclP->GL,lclP->GL_U,lclP->GL_V);
488:       LCLScatter(lclP,lclP->GAugL,lclP->GAugL_U,lclP->GAugL_V);
489:       TaoComputeObjectiveAndGradient(tao,tao->solution,&f,tao->gradient);
490:       LCLScatter(lclP,tao->gradient,lclP->GU,lclP->GV);

492:       /* Compute the reduced gradient at the new point */

494:       TaoComputeJacobianState(tao,lclP->X0,tao->jacobian_state,tao->jacobian_state_pre,tao->jacobian_state_inv);
495:       TaoComputeJacobianDesign(tao,lclP->X0,tao->jacobian_design);

497:       /* p2 */
498:       /* Compute multipliers, using lamda-rho*con as an initial guess in PCG */
499:       if (phase2_iter==0){
500:         VecWAXPY(lclP->lamda,-lclP->rho,lclP->con1,lclP->lamda0);
501:       } else {
502:         VecAXPY(lclP->lamda,-lclP->rho,tao->constraints);
503:       }

505:       MatIsSymmetricKnown(tao->jacobian_state,&set,&flag);
506:       if (tao->jacobian_state_pre) {
507:         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;

514:       lclP->solve_type = LCL_ADJOINT2;
515:       if (tao->jacobian_state_inv) {
516:         if (symmetric) {
517:           MatMult(tao->jacobian_state_inv, lclP->GU, lclP->lamda);
518:         } else {
519:           MatMultTranspose(tao->jacobian_state_inv, lclP->GU, lclP->lamda);
520:         }
521:       } else {
522:         if (symmetric) {
523:           KSPSolve(tao->ksp, lclP->GU,  lclP->lamda);
524:         } else {
525:           KSPSolveTranspose(tao->ksp, lclP->GU,  lclP->lamda);
526:         }
527:         KSPGetIterationNumber(tao->ksp,&its);
528:         tao->ksp_its += its;
529:         tao->ksp_tot_its += its;
530:       }

532:       MatMultTranspose(tao->jacobian_design,lclP->lamda,lclP->g2);
533:       VecAXPY(lclP->g2,-1.0,lclP->GV);

535:       VecScale(lclP->g2,-1.0);

537:       /* Update the quasi-newton approximation */
538:       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 */

541:     }

543:     VecCopy(lclP->lamda,lclP->lamda0);

545:     /* Evaluate Function, Gradient, Constraints, and Jacobian */
546:     TaoComputeObjectiveAndGradient(tao,tao->solution,&f,tao->gradient);
547:     LCLScatter(lclP,tao->solution,lclP->U,lclP->V);
548:     LCLScatter(lclP,tao->gradient,lclP->GU,lclP->GV);

550:     TaoComputeJacobianState(tao,tao->solution,tao->jacobian_state,tao->jacobian_state_pre,tao->jacobian_state_inv);
551:     TaoComputeJacobianDesign(tao,tao->solution,tao->jacobian_design);
552:     TaoComputeConstraints(tao,tao->solution, tao->constraints);

554:     LCLComputeAugmentedLagrangianAndGradient(tao->linesearch,tao->solution,&lclP->aug,lclP->GAugL,tao);

556:     VecNorm(lclP->GAugL, NORM_2, &mnorm);

558:     /* Evaluate constraint norm */
559:     VecNorm(tao->constraints, NORM_2, &cnorm);

561:     /* Monitor convergence */
562:     tao->niter++;
563:     TaoLogConvergenceHistory(tao,f,mnorm,cnorm,tao->ksp_its);
564:     TaoMonitor(tao,tao->niter,f,mnorm,cnorm,step);
565:     (*tao->ops->convergencetest)(tao,tao->cnvP);
566:   }
567:   return(0);
568: }

570: /*MC
571:  TAOLCL - linearly constrained lagrangian method for pde-constrained optimization

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

584:   Level: beginner
585: M*/
586: PETSC_EXTERN PetscErrorCode TaoCreate_LCL(Tao tao)
587: {
588:   TAO_LCL        *lclP;
590:   const char     *morethuente_type = TAOLINESEARCHMT;

593:   tao->ops->setup = TaoSetup_LCL;
594:   tao->ops->solve = TaoSolve_LCL;
595:   tao->ops->view = TaoView_LCL;
596:   tao->ops->setfromoptions = TaoSetFromOptions_LCL;
597:   tao->ops->destroy = TaoDestroy_LCL;
598:   PetscNewLog(tao,&lclP);
599:   tao->data = (void*)lclP;

601:   /* Override default settings (unless already changed) */
602:   if (!tao->max_it_changed) tao->max_it = 200;
603:   if (!tao->catol_changed) tao->catol = 1.0e-4;
604:   if (!tao->gatol_changed) tao->gttol = 1.0e-4;
605:   if (!tao->grtol_changed) tao->gttol = 1.0e-4;
606:   if (!tao->gttol_changed) tao->gttol = 1.0e-4;
607:   lclP->rho0 = 1.0e-4;
608:   lclP->rhomax=1e5;
609:   lclP->eps1 = 1.0e-8;
610:   lclP->eps2 = 0.0;
611:   lclP->solve_type=2;
612:   lclP->tau[0] = lclP->tau[1] = lclP->tau[2] = lclP->tau[3] = 1.0e-4;
613:   TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch);
614:   PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1);
615:   TaoLineSearchSetType(tao->linesearch, morethuente_type);
616:   TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix);

618:   TaoLineSearchSetObjectiveAndGradientRoutine(tao->linesearch,LCLComputeAugmentedLagrangianAndGradient, tao);
619:   KSPCreate(((PetscObject)tao)->comm,&tao->ksp);
620:   PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1);
621:   KSPSetOptionsPrefix(tao->ksp, tao->hdr.prefix);
622:   KSPSetFromOptions(tao->ksp);
623: 
624:   MatCreate(((PetscObject)tao)->comm, &lclP->R);
625:   MatSetType(lclP->R, MATLMVMBFGS);
626:   return(0);
627: }

629: static PetscErrorCode LCLComputeLagrangianAndGradient(TaoLineSearch ls, Vec X, PetscReal *f, Vec G, void *ptr)
630: {
631:   Tao            tao = (Tao)ptr;
632:   TAO_LCL        *lclP = (TAO_LCL*)tao->data;
633:   PetscBool      set,pset,flag,pflag,symmetric;
634:   PetscReal      cdotl;

638:   TaoComputeObjectiveAndGradient(tao,X,f,G);
639:   LCLScatter(lclP,G,lclP->GU,lclP->GV);
640:   if (lclP->recompute_jacobian_flag) {
641:     TaoComputeJacobianState(tao,X,tao->jacobian_state,tao->jacobian_state_pre,tao->jacobian_state_inv);
642:     TaoComputeJacobianDesign(tao,X,tao->jacobian_design);
643:   }
644:   TaoComputeConstraints(tao,X, tao->constraints);
645:   MatIsSymmetricKnown(tao->jacobian_state,&set,&flag);
646:   if (tao->jacobian_state_pre) {
647:     MatIsSymmetricKnown(tao->jacobian_state_pre,&pset,&pflag);
648:   } else {
649:     pset = pflag = PETSC_TRUE;
650:   }
651:   if (set && pset && flag && pflag) symmetric = PETSC_TRUE;
652:   else symmetric = PETSC_FALSE;

654:   VecDot(lclP->lamda0, tao->constraints, &cdotl);
655:   lclP->lgn = *f - cdotl;

657:   /* Gradient of Lagrangian GL = G - J' * lamda */
658:   /*      WU = A' * WL
659:           WV = B' * WL */
660:   if (symmetric) {
661:     MatMult(tao->jacobian_state,lclP->lamda0,lclP->GL_U);
662:   } else {
663:     MatMultTranspose(tao->jacobian_state,lclP->lamda0,lclP->GL_U);
664:   }
665:   MatMultTranspose(tao->jacobian_design,lclP->lamda0,lclP->GL_V);
666:   VecScale(lclP->GL_U,-1.0);
667:   VecScale(lclP->GL_V,-1.0);
668:   VecAXPY(lclP->GL_U,1.0,lclP->GU);
669:   VecAXPY(lclP->GL_V,1.0,lclP->GV);
670:   LCLGather(lclP,lclP->GL_U,lclP->GL_V,lclP->GL);

672:   f[0] = lclP->lgn;
673:   VecCopy(lclP->GL,G);
674:   return(0);
675: }

677: static PetscErrorCode LCLComputeAugmentedLagrangianAndGradient(TaoLineSearch ls, Vec X, PetscReal *f, Vec G, void *ptr)
678: {
679:   Tao            tao = (Tao)ptr;
680:   TAO_LCL        *lclP = (TAO_LCL*)tao->data;
681:   PetscReal      con2;
682:   PetscBool      flag,pflag,set,pset,symmetric;

686:   LCLComputeLagrangianAndGradient(tao->linesearch,X,f,G,tao);
687:   LCLScatter(lclP,G,lclP->GL_U,lclP->GL_V);
688:   VecDot(tao->constraints,tao->constraints,&con2);
689:   lclP->aug = lclP->lgn + 0.5*lclP->rho*con2;

691:   /* Gradient of Aug. Lagrangian GAugL = GL + rho * J' c */
692:   /*      WU = A' * c
693:           WV = B' * c */
694:   MatIsSymmetricKnown(tao->jacobian_state,&set,&flag);
695:   if (tao->jacobian_state_pre) {
696:     MatIsSymmetricKnown(tao->jacobian_state_pre,&pset,&pflag);
697:   } else {
698:     pset = pflag = PETSC_TRUE;
699:   }
700:   if (set && pset && flag && pflag) symmetric = PETSC_TRUE;
701:   else symmetric = PETSC_FALSE;

703:   if (symmetric) {
704:     MatMult(tao->jacobian_state,tao->constraints,lclP->GAugL_U);
705:   } else {
706:     MatMultTranspose(tao->jacobian_state,tao->constraints,lclP->GAugL_U);
707:   }

709:   MatMultTranspose(tao->jacobian_design,tao->constraints,lclP->GAugL_V);
710:   VecAYPX(lclP->GAugL_U,lclP->rho,lclP->GL_U);
711:   VecAYPX(lclP->GAugL_V,lclP->rho,lclP->GL_V);
712:   LCLGather(lclP,lclP->GAugL_U,lclP->GAugL_V,lclP->GAugL);

714:   f[0] = lclP->aug;
715:   VecCopy(lclP->GAugL,G);
716:   return(0);
717: }

719: PetscErrorCode LCLGather(TAO_LCL *lclP, Vec u, Vec v, Vec x)
720: {
723:   VecScatterBegin(lclP->state_scatter, u, x, INSERT_VALUES, SCATTER_REVERSE);
724:   VecScatterEnd(lclP->state_scatter, u, x, INSERT_VALUES, SCATTER_REVERSE);
725:   VecScatterBegin(lclP->design_scatter, v, x, INSERT_VALUES, SCATTER_REVERSE);
726:   VecScatterEnd(lclP->design_scatter, v, x, INSERT_VALUES, SCATTER_REVERSE);
727:   return(0);

729: }
730: PetscErrorCode LCLScatter(TAO_LCL *lclP, Vec x, Vec u, Vec v)
731: {
734:   VecScatterBegin(lclP->state_scatter, x, u, INSERT_VALUES, SCATTER_FORWARD);
735:   VecScatterEnd(lclP->state_scatter, x, u, INSERT_VALUES, SCATTER_FORWARD);
736:   VecScatterBegin(lclP->design_scatter, x, v, INSERT_VALUES, SCATTER_FORWARD);
737:   VecScatterEnd(lclP->design_scatter, x, v, INSERT_VALUES, SCATTER_FORWARD);
738:   return(0);

740: }