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