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) PetscCall((*tao->ops->update)(tao, tao->niter, tao->user_update)); 250 tao->ksp_its=0; 251 /* Compute a descent direction for the linearly constrained subproblem 252 minimize f(u+du, v+dv) 253 s.t. A(u0,v0)du + B(u0,v0)dv = -g(u0,v0) */ 254 255 /* Store the points around the linearization */ 256 PetscCall(VecCopy(lclP->U, lclP->U0)); 257 PetscCall(VecCopy(lclP->V, lclP->V0)); 258 PetscCall(VecCopy(lclP->GU,lclP->GU0)); 259 PetscCall(VecCopy(lclP->GV,lclP->GV0)); 260 PetscCall(VecCopy(lclP->GAugL_U,lclP->GAugL_U0)); 261 PetscCall(VecCopy(lclP->GAugL_V,lclP->GAugL_V0)); 262 PetscCall(VecCopy(lclP->GL_U,lclP->GL_U0)); 263 PetscCall(VecCopy(lclP->GL_V,lclP->GL_V0)); 264 265 lclP->aug0 = lclP->aug; 266 lclP->lgn0 = lclP->lgn; 267 268 /* Given the design variables, we need to project the current iterate 269 onto the linearized constraint. We choose to fix the design variables 270 and solve the linear system for the state variables. The resulting 271 point is the Newton direction */ 272 273 /* Solve r = A\con */ 274 lclP->solve_type = LCL_FORWARD1; 275 PetscCall(VecSet(lclP->r,0.0)); /* Initial guess in CG */ 276 277 if (tao->jacobian_state_inv) { 278 PetscCall(MatMult(tao->jacobian_state_inv, tao->constraints, lclP->r)); 279 } else { 280 PetscCall(KSPSetOperators(tao->ksp, tao->jacobian_state, tao->jacobian_state_pre)); 281 PetscCall(KSPSolve(tao->ksp, tao->constraints, lclP->r)); 282 PetscCall(KSPGetIterationNumber(tao->ksp,&its)); 283 tao->ksp_its+=its; 284 tao->ksp_tot_its+=tao->ksp_its; 285 } 286 287 /* Set design step direction dv to zero */ 288 PetscCall(VecSet(lclP->s, 0.0)); 289 290 /* 291 Check sufficient descent for constraint merit function .5*||con||^2 292 con' Ak r >= eps1 ||r||^(2+eps2) 293 */ 294 295 /* Compute WU= Ak' * con */ 296 if (symmetric) { 297 PetscCall(MatMult(tao->jacobian_state,tao->constraints,lclP->WU)); 298 } else { 299 PetscCall(MatMultTranspose(tao->jacobian_state,tao->constraints,lclP->WU)); 300 } 301 /* Compute r * Ak' * con */ 302 PetscCall(VecDot(lclP->r,lclP->WU,&rWU)); 303 304 /* compute ||r||^(2+eps2) */ 305 PetscCall(VecNorm(lclP->r,NORM_2,&r2)); 306 r2 = PetscPowScalar(r2,2.0+lclP->eps2); 307 adec = lclP->eps1 * r2; 308 309 if (rWU < adec) { 310 PetscCall(PetscInfo(tao,"Newton direction not descent for constraint, feasibility phase required\n")); 311 if (lclP->verbose) { 312 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Newton direction not descent for constraint: %g -- using steepest descent\n",(double)descent)); 313 } 314 315 PetscCall(PetscInfo(tao,"Using steepest descent direction instead.\n")); 316 PetscCall(VecSet(lclP->r,0.0)); 317 PetscCall(VecAXPY(lclP->r,-1.0,lclP->WU)); 318 PetscCall(VecDot(lclP->r,lclP->r,&rWU)); 319 PetscCall(VecNorm(lclP->r,NORM_2,&r2)); 320 r2 = PetscPowScalar(r2,2.0+lclP->eps2); 321 PetscCall(VecDot(lclP->r,lclP->GAugL_U,&descent)); 322 adec = lclP->eps1 * r2; 323 } 324 325 /* 326 Check descent for aug. lagrangian 327 r' (GUk - Ak'*yk - rho*Ak'*con) <= -eps1 ||r||^(2+eps2) 328 GL_U = GUk - Ak'*yk 329 WU = Ak'*con 330 adec=eps1||r||^(2+eps2) 331 332 ==> 333 Check r'GL_U - rho*r'WU <= adec 334 */ 335 336 PetscCall(VecDot(lclP->r,lclP->GL_U,&rGL_U)); 337 aldescent = rGL_U - lclP->rho*rWU; 338 if (aldescent > -adec) { 339 if (lclP->verbose) { 340 PetscCall(PetscPrintf(PETSC_COMM_WORLD," Newton direction not descent for augmented Lagrangian: %g",(double)aldescent)); 341 } 342 PetscCall(PetscInfo(tao,"Newton direction not descent for augmented Lagrangian: %g\n",(double)aldescent)); 343 lclP->rho = (rGL_U - adec)/rWU; 344 if (lclP->rho > lclP->rhomax) { 345 lclP->rho = lclP->rhomax; 346 SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"rho=%g > rhomax, case not implemented. Increase rhomax (-tao_lcl_rhomax)",(double)lclP->rho); 347 } 348 if (lclP->verbose) { 349 PetscCall(PetscPrintf(PETSC_COMM_WORLD," Increasing penalty parameter to %g\n",(double)lclP->rho)); 350 } 351 PetscCall(PetscInfo(tao," Increasing penalty parameter to %g\n",(double)lclP->rho)); 352 } 353 354 PetscCall(LCLComputeAugmentedLagrangianAndGradient(tao->linesearch,tao->solution,&lclP->aug,lclP->GAugL,tao)); 355 PetscCall(LCLScatter(lclP,lclP->GL,lclP->GL_U,lclP->GL_V)); 356 PetscCall(LCLScatter(lclP,lclP->GAugL,lclP->GAugL_U,lclP->GAugL_V)); 357 358 /* We now minimize the augmented Lagrangian along the Newton direction */ 359 PetscCall(VecScale(lclP->r,-1.0)); 360 PetscCall(LCLGather(lclP, lclP->r,lclP->s,tao->stepdirection)); 361 PetscCall(VecScale(lclP->r,-1.0)); 362 PetscCall(LCLGather(lclP, lclP->GAugL_U0, lclP->GAugL_V0, lclP->GAugL)); 363 PetscCall(LCLGather(lclP, lclP->U0,lclP->V0,lclP->X0)); 364 365 lclP->recompute_jacobian_flag = PETSC_TRUE; 366 367 PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch,1.0)); 368 PetscCall(TaoLineSearchSetType(tao->linesearch, TAOLINESEARCHMT)); 369 PetscCall(TaoLineSearchSetFromOptions(tao->linesearch)); 370 PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &lclP->aug, lclP->GAugL, tao->stepdirection, &step, &ls_reason)); 371 if (lclP->verbose) { 372 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Steplength = %g\n",(double)step)); 373 } 374 375 PetscCall(LCLScatter(lclP,tao->solution,lclP->U,lclP->V)); 376 PetscCall(TaoComputeObjectiveAndGradient(tao,tao->solution,&f,tao->gradient)); 377 PetscCall(LCLScatter(lclP,tao->gradient,lclP->GU,lclP->GV)); 378 379 PetscCall(LCLScatter(lclP,lclP->GAugL,lclP->GAugL_U,lclP->GAugL_V)); 380 381 /* Check convergence */ 382 PetscCall(VecNorm(lclP->GAugL, NORM_2, &mnorm)); 383 PetscCall(VecNorm(tao->constraints, NORM_2, &cnorm)); 384 PetscCall(TaoLogConvergenceHistory(tao,f,mnorm,cnorm,tao->ksp_its)); 385 PetscCall(TaoMonitor(tao,tao->niter,f,mnorm,cnorm,step)); 386 PetscCall((*tao->ops->convergencetest)(tao,tao->cnvP)); 387 if (tao->reason != TAO_CONTINUE_ITERATING) { 388 break; 389 } 390 391 /* TODO: use a heuristic to choose how many iterations should be performed within phase 2 */ 392 for (phase2_iter=0; phase2_iter<lclP->phase2_niter; phase2_iter++) { 393 /* We now minimize the objective function starting from the fraction of 394 the Newton point accepted by applying one step of a reduced-space 395 method. The optimization problem is: 396 397 minimize f(u+du, v+dv) 398 s. t. A(u0,v0)du + B(u0,v0)du = -alpha g(u0,v0) 399 400 In particular, we have that 401 du = -inv(A)*(Bdv + alpha g) */ 402 403 PetscCall(TaoComputeJacobianState(tao,lclP->X0,tao->jacobian_state,tao->jacobian_state_pre,tao->jacobian_state_inv)); 404 PetscCall(TaoComputeJacobianDesign(tao,lclP->X0,tao->jacobian_design)); 405 406 /* Store V and constraints */ 407 PetscCall(VecCopy(lclP->V, lclP->V1)); 408 PetscCall(VecCopy(tao->constraints,lclP->con1)); 409 410 /* Compute multipliers */ 411 /* p1 */ 412 PetscCall(VecSet(lclP->lamda,0.0)); /* Initial guess in CG */ 413 lclP->solve_type = LCL_ADJOINT1; 414 PetscCall(MatIsSymmetricKnown(tao->jacobian_state,&set,&flag)); 415 if (tao->jacobian_state_pre) { 416 PetscCall(MatIsSymmetricKnown(tao->jacobian_state_pre,&pset,&pflag)); 417 } else { 418 pset = pflag = PETSC_TRUE; 419 } 420 if (set && pset && flag && pflag) symmetric = PETSC_TRUE; 421 else symmetric = PETSC_FALSE; 422 423 if (tao->jacobian_state_inv) { 424 if (symmetric) { 425 PetscCall(MatMult(tao->jacobian_state_inv, lclP->GAugL_U, lclP->lamda)); 426 } else { 427 PetscCall(MatMultTranspose(tao->jacobian_state_inv, lclP->GAugL_U, lclP->lamda)); 428 } 429 } else { 430 if (symmetric) { 431 PetscCall(KSPSolve(tao->ksp, lclP->GAugL_U, lclP->lamda)); 432 } else { 433 PetscCall(KSPSolveTranspose(tao->ksp, lclP->GAugL_U, lclP->lamda)); 434 } 435 PetscCall(KSPGetIterationNumber(tao->ksp,&its)); 436 tao->ksp_its+=its; 437 tao->ksp_tot_its+=its; 438 } 439 PetscCall(MatMultTranspose(tao->jacobian_design,lclP->lamda,lclP->g1)); 440 PetscCall(VecAXPY(lclP->g1,-1.0,lclP->GAugL_V)); 441 442 /* Compute the limited-memory quasi-newton direction */ 443 if (tao->niter > 0) { 444 PetscCall(MatSolve(lclP->R,lclP->g1,lclP->s)); 445 PetscCall(VecDot(lclP->s,lclP->g1,&descent)); 446 if (descent <= 0) { 447 if (lclP->verbose) { 448 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Reduced-space direction not descent: %g\n",(double)descent)); 449 } 450 PetscCall(VecCopy(lclP->g1,lclP->s)); 451 } 452 } else { 453 PetscCall(VecCopy(lclP->g1,lclP->s)); 454 } 455 PetscCall(VecScale(lclP->g1,-1.0)); 456 457 /* Recover the full space direction */ 458 PetscCall(MatMult(tao->jacobian_design,lclP->s,lclP->WU)); 459 /* PetscCall(VecSet(lclP->r,0.0)); */ /* Initial guess in CG */ 460 lclP->solve_type = LCL_FORWARD2; 461 if (tao->jacobian_state_inv) { 462 PetscCall(MatMult(tao->jacobian_state_inv,lclP->WU,lclP->r)); 463 } else { 464 PetscCall(KSPSolve(tao->ksp, lclP->WU, lclP->r)); 465 PetscCall(KSPGetIterationNumber(tao->ksp,&its)); 466 tao->ksp_its += its; 467 tao->ksp_tot_its+=its; 468 } 469 470 /* We now minimize the augmented Lagrangian along the direction -r,s */ 471 PetscCall(VecScale(lclP->r, -1.0)); 472 PetscCall(LCLGather(lclP,lclP->r,lclP->s,tao->stepdirection)); 473 PetscCall(VecScale(lclP->r, -1.0)); 474 lclP->recompute_jacobian_flag = PETSC_TRUE; 475 476 PetscCall(TaoLineSearchSetInitialStepLength(tao->linesearch,1.0)); 477 PetscCall(TaoLineSearchSetType(tao->linesearch,TAOLINESEARCHMT)); 478 PetscCall(TaoLineSearchSetFromOptions(tao->linesearch)); 479 PetscCall(TaoLineSearchApply(tao->linesearch, tao->solution, &lclP->aug, lclP->GAugL, tao->stepdirection,&step,&ls_reason)); 480 if (lclP->verbose) { 481 PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Reduced-space steplength = %g\n",(double)step)); 482 } 483 484 PetscCall(LCLScatter(lclP,tao->solution,lclP->U,lclP->V)); 485 PetscCall(LCLScatter(lclP,lclP->GL,lclP->GL_U,lclP->GL_V)); 486 PetscCall(LCLScatter(lclP,lclP->GAugL,lclP->GAugL_U,lclP->GAugL_V)); 487 PetscCall(TaoComputeObjectiveAndGradient(tao,tao->solution,&f,tao->gradient)); 488 PetscCall(LCLScatter(lclP,tao->gradient,lclP->GU,lclP->GV)); 489 490 /* Compute the reduced gradient at the new point */ 491 492 PetscCall(TaoComputeJacobianState(tao,lclP->X0,tao->jacobian_state,tao->jacobian_state_pre,tao->jacobian_state_inv)); 493 PetscCall(TaoComputeJacobianDesign(tao,lclP->X0,tao->jacobian_design)); 494 495 /* p2 */ 496 /* Compute multipliers, using lamda-rho*con as an initial guess in PCG */ 497 if (phase2_iter==0) { 498 PetscCall(VecWAXPY(lclP->lamda,-lclP->rho,lclP->con1,lclP->lamda0)); 499 } else { 500 PetscCall(VecAXPY(lclP->lamda,-lclP->rho,tao->constraints)); 501 } 502 503 PetscCall(MatIsSymmetricKnown(tao->jacobian_state,&set,&flag)); 504 if (tao->jacobian_state_pre) { 505 PetscCall(MatIsSymmetricKnown(tao->jacobian_state_pre,&pset,&pflag)); 506 } else { 507 pset = pflag = PETSC_TRUE; 508 } 509 if (set && pset && flag && pflag) symmetric = PETSC_TRUE; 510 else symmetric = PETSC_FALSE; 511 512 lclP->solve_type = LCL_ADJOINT2; 513 if (tao->jacobian_state_inv) { 514 if (symmetric) { 515 PetscCall(MatMult(tao->jacobian_state_inv, lclP->GU, lclP->lamda)); 516 } else { 517 PetscCall(MatMultTranspose(tao->jacobian_state_inv, lclP->GU, lclP->lamda)); 518 } 519 } else { 520 if (symmetric) { 521 PetscCall(KSPSolve(tao->ksp, lclP->GU, lclP->lamda)); 522 } else { 523 PetscCall(KSPSolveTranspose(tao->ksp, lclP->GU, lclP->lamda)); 524 } 525 PetscCall(KSPGetIterationNumber(tao->ksp,&its)); 526 tao->ksp_its += its; 527 tao->ksp_tot_its += its; 528 } 529 530 PetscCall(MatMultTranspose(tao->jacobian_design,lclP->lamda,lclP->g2)); 531 PetscCall(VecAXPY(lclP->g2,-1.0,lclP->GV)); 532 533 PetscCall(VecScale(lclP->g2,-1.0)); 534 535 /* Update the quasi-newton approximation */ 536 PetscCall(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 */ 538 539 } 540 541 PetscCall(VecCopy(lclP->lamda,lclP->lamda0)); 542 543 /* Evaluate Function, Gradient, Constraints, and Jacobian */ 544 PetscCall(TaoComputeObjectiveAndGradient(tao,tao->solution,&f,tao->gradient)); 545 PetscCall(LCLScatter(lclP,tao->solution,lclP->U,lclP->V)); 546 PetscCall(LCLScatter(lclP,tao->gradient,lclP->GU,lclP->GV)); 547 548 PetscCall(TaoComputeJacobianState(tao,tao->solution,tao->jacobian_state,tao->jacobian_state_pre,tao->jacobian_state_inv)); 549 PetscCall(TaoComputeJacobianDesign(tao,tao->solution,tao->jacobian_design)); 550 PetscCall(TaoComputeConstraints(tao,tao->solution, tao->constraints)); 551 552 PetscCall(LCLComputeAugmentedLagrangianAndGradient(tao->linesearch,tao->solution,&lclP->aug,lclP->GAugL,tao)); 553 554 PetscCall(VecNorm(lclP->GAugL, NORM_2, &mnorm)); 555 556 /* Evaluate constraint norm */ 557 PetscCall(VecNorm(tao->constraints, NORM_2, &cnorm)); 558 559 /* Monitor convergence */ 560 tao->niter++; 561 PetscCall(TaoLogConvergenceHistory(tao,f,mnorm,cnorm,tao->ksp_its)); 562 PetscCall(TaoMonitor(tao,tao->niter,f,mnorm,cnorm,step)); 563 PetscCall((*tao->ops->convergencetest)(tao,tao->cnvP)); 564 } 565 PetscFunctionReturn(0); 566 } 567 568 /*MC 569 TAOLCL - linearly constrained lagrangian method for pde-constrained optimization 570 571 + -tao_lcl_eps1 - epsilon 1 tolerance 572 . -tao_lcl_eps2","epsilon 2 tolerance","",lclP->eps2,&lclP->eps2,NULL); 573 . -tao_lcl_rho0","init value for rho","",lclP->rho0,&lclP->rho0,NULL); 574 . -tao_lcl_rhomax","max value for rho","",lclP->rhomax,&lclP->rhomax,NULL); 575 . -tao_lcl_phase2_niter - Number of phase 2 iterations in LCL algorithm 576 . -tao_lcl_verbose - Print verbose output if True 577 . -tao_lcl_tola - Tolerance for first forward solve 578 . -tao_lcl_tolb - Tolerance for first adjoint solve 579 . -tao_lcl_tolc - Tolerance for second forward solve 580 - -tao_lcl_told - Tolerance for second adjoint solve 581 582 Level: beginner 583 M*/ 584 PETSC_EXTERN PetscErrorCode TaoCreate_LCL(Tao tao) 585 { 586 TAO_LCL *lclP; 587 const char *morethuente_type = TAOLINESEARCHMT; 588 589 PetscFunctionBegin; 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 PetscCall(PetscNewLog(tao,&lclP)); 596 tao->data = (void*)lclP; 597 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 PetscCall(TaoLineSearchCreate(((PetscObject)tao)->comm, &tao->linesearch)); 611 PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->linesearch, (PetscObject)tao, 1)); 612 PetscCall(TaoLineSearchSetType(tao->linesearch, morethuente_type)); 613 PetscCall(TaoLineSearchSetOptionsPrefix(tao->linesearch,tao->hdr.prefix)); 614 615 PetscCall(TaoLineSearchSetObjectiveAndGradientRoutine(tao->linesearch,LCLComputeAugmentedLagrangianAndGradient, tao)); 616 PetscCall(KSPCreate(((PetscObject)tao)->comm,&tao->ksp)); 617 PetscCall(PetscObjectIncrementTabLevel((PetscObject)tao->ksp, (PetscObject)tao, 1)); 618 PetscCall(KSPSetOptionsPrefix(tao->ksp, tao->hdr.prefix)); 619 PetscCall(KSPSetFromOptions(tao->ksp)); 620 621 PetscCall(MatCreate(((PetscObject)tao)->comm, &lclP->R)); 622 PetscCall(MatSetType(lclP->R, MATLMVMBFGS)); 623 PetscFunctionReturn(0); 624 } 625 626 static PetscErrorCode LCLComputeLagrangianAndGradient(TaoLineSearch ls, Vec X, PetscReal *f, Vec G, void *ptr) 627 { 628 Tao tao = (Tao)ptr; 629 TAO_LCL *lclP = (TAO_LCL*)tao->data; 630 PetscBool set,pset,flag,pflag,symmetric; 631 PetscReal cdotl; 632 633 PetscFunctionBegin; 634 PetscCall(TaoComputeObjectiveAndGradient(tao,X,f,G)); 635 PetscCall(LCLScatter(lclP,G,lclP->GU,lclP->GV)); 636 if (lclP->recompute_jacobian_flag) { 637 PetscCall(TaoComputeJacobianState(tao,X,tao->jacobian_state,tao->jacobian_state_pre,tao->jacobian_state_inv)); 638 PetscCall(TaoComputeJacobianDesign(tao,X,tao->jacobian_design)); 639 } 640 PetscCall(TaoComputeConstraints(tao,X, tao->constraints)); 641 PetscCall(MatIsSymmetricKnown(tao->jacobian_state,&set,&flag)); 642 if (tao->jacobian_state_pre) { 643 PetscCall(MatIsSymmetricKnown(tao->jacobian_state_pre,&pset,&pflag)); 644 } else { 645 pset = pflag = PETSC_TRUE; 646 } 647 if (set && pset && flag && pflag) symmetric = PETSC_TRUE; 648 else symmetric = PETSC_FALSE; 649 650 PetscCall(VecDot(lclP->lamda0, tao->constraints, &cdotl)); 651 lclP->lgn = *f - cdotl; 652 653 /* Gradient of Lagrangian GL = G - J' * lamda */ 654 /* WU = A' * WL 655 WV = B' * WL */ 656 if (symmetric) { 657 PetscCall(MatMult(tao->jacobian_state,lclP->lamda0,lclP->GL_U)); 658 } else { 659 PetscCall(MatMultTranspose(tao->jacobian_state,lclP->lamda0,lclP->GL_U)); 660 } 661 PetscCall(MatMultTranspose(tao->jacobian_design,lclP->lamda0,lclP->GL_V)); 662 PetscCall(VecScale(lclP->GL_U,-1.0)); 663 PetscCall(VecScale(lclP->GL_V,-1.0)); 664 PetscCall(VecAXPY(lclP->GL_U,1.0,lclP->GU)); 665 PetscCall(VecAXPY(lclP->GL_V,1.0,lclP->GV)); 666 PetscCall(LCLGather(lclP,lclP->GL_U,lclP->GL_V,lclP->GL)); 667 668 f[0] = lclP->lgn; 669 PetscCall(VecCopy(lclP->GL,G)); 670 PetscFunctionReturn(0); 671 } 672 673 static PetscErrorCode LCLComputeAugmentedLagrangianAndGradient(TaoLineSearch ls, Vec X, PetscReal *f, Vec G, void *ptr) 674 { 675 Tao tao = (Tao)ptr; 676 TAO_LCL *lclP = (TAO_LCL*)tao->data; 677 PetscReal con2; 678 PetscBool flag,pflag,set,pset,symmetric; 679 680 PetscFunctionBegin; 681 PetscCall(LCLComputeLagrangianAndGradient(tao->linesearch,X,f,G,tao)); 682 PetscCall(LCLScatter(lclP,G,lclP->GL_U,lclP->GL_V)); 683 PetscCall(VecDot(tao->constraints,tao->constraints,&con2)); 684 lclP->aug = lclP->lgn + 0.5*lclP->rho*con2; 685 686 /* Gradient of Aug. Lagrangian GAugL = GL + rho * J' c */ 687 /* WU = A' * c 688 WV = B' * c */ 689 PetscCall(MatIsSymmetricKnown(tao->jacobian_state,&set,&flag)); 690 if (tao->jacobian_state_pre) { 691 PetscCall(MatIsSymmetricKnown(tao->jacobian_state_pre,&pset,&pflag)); 692 } else { 693 pset = pflag = PETSC_TRUE; 694 } 695 if (set && pset && flag && pflag) symmetric = PETSC_TRUE; 696 else symmetric = PETSC_FALSE; 697 698 if (symmetric) { 699 PetscCall(MatMult(tao->jacobian_state,tao->constraints,lclP->GAugL_U)); 700 } else { 701 PetscCall(MatMultTranspose(tao->jacobian_state,tao->constraints,lclP->GAugL_U)); 702 } 703 704 PetscCall(MatMultTranspose(tao->jacobian_design,tao->constraints,lclP->GAugL_V)); 705 PetscCall(VecAYPX(lclP->GAugL_U,lclP->rho,lclP->GL_U)); 706 PetscCall(VecAYPX(lclP->GAugL_V,lclP->rho,lclP->GL_V)); 707 PetscCall(LCLGather(lclP,lclP->GAugL_U,lclP->GAugL_V,lclP->GAugL)); 708 709 f[0] = lclP->aug; 710 PetscCall(VecCopy(lclP->GAugL,G)); 711 PetscFunctionReturn(0); 712 } 713 714 PetscErrorCode LCLGather(TAO_LCL *lclP, Vec u, Vec v, Vec x) 715 { 716 PetscFunctionBegin; 717 PetscCall(VecScatterBegin(lclP->state_scatter, u, x, INSERT_VALUES, SCATTER_REVERSE)); 718 PetscCall(VecScatterEnd(lclP->state_scatter, u, x, INSERT_VALUES, SCATTER_REVERSE)); 719 PetscCall(VecScatterBegin(lclP->design_scatter, v, x, INSERT_VALUES, SCATTER_REVERSE)); 720 PetscCall(VecScatterEnd(lclP->design_scatter, v, x, INSERT_VALUES, SCATTER_REVERSE)); 721 PetscFunctionReturn(0); 722 723 } 724 PetscErrorCode LCLScatter(TAO_LCL *lclP, Vec x, Vec u, Vec v) 725 { 726 PetscFunctionBegin; 727 PetscCall(VecScatterBegin(lclP->state_scatter, x, u, INSERT_VALUES, SCATTER_FORWARD)); 728 PetscCall(VecScatterEnd(lclP->state_scatter, x, u, INSERT_VALUES, SCATTER_FORWARD)); 729 PetscCall(VecScatterBegin(lclP->design_scatter, x, v, INSERT_VALUES, SCATTER_FORWARD)); 730 PetscCall(VecScatterEnd(lclP->design_scatter, x, v, INSERT_VALUES, SCATTER_FORWARD)); 731 PetscFunctionReturn(0); 732 733 } 734