1 #include <../src/tao/constrained/impls/almm/almm.h> /*I "petsctao.h" I*/ /*I "petscvec.h" I*/ 2 #include <petsctao.h> 3 #include <petsc/private/petscimpl.h> 4 #include <petsc/private/vecimpl.h> 5 6 static PetscErrorCode TaoALMMCombinePrimal_Private(Tao,Vec,Vec,Vec); 7 static PetscErrorCode TaoALMMCombineDual_Private(Tao,Vec,Vec,Vec); 8 static PetscErrorCode TaoALMMSplitPrimal_Private(Tao,Vec,Vec,Vec); 9 static PetscErrorCode TaoALMMComputeOptimalityNorms_Private(Tao); 10 static PetscErrorCode TaoALMMComputeAugLagAndGradient_Private(Tao); 11 static PetscErrorCode TaoALMMComputePHRLagAndGradient_Private(Tao); 12 13 static PetscErrorCode TaoSolve_ALMM(Tao tao) 14 { 15 TAO_ALMM *auglag = (TAO_ALMM*)tao->data; 16 TaoConvergedReason reason; 17 PetscReal updated; 18 19 PetscFunctionBegin; 20 /* reset initial multiplier/slack guess */ 21 if (!tao->recycle) { 22 if (tao->ineq_constrained) { 23 PetscCall(VecZeroEntries(auglag->Ps)); 24 PetscCall(TaoALMMCombinePrimal_Private(tao, auglag->Px, auglag->Ps, auglag->P)); 25 PetscCall(VecZeroEntries(auglag->Yi)); 26 } 27 if (tao->eq_constrained) { 28 PetscCall(VecZeroEntries(auglag->Ye)); 29 } 30 } 31 32 /* compute initial nonlinear Lagrangian and its derivatives */ 33 PetscCall((*auglag->sub_obj)(tao)); 34 PetscCall(TaoALMMComputeOptimalityNorms_Private(tao)); 35 /* print initial step and check convergence */ 36 PetscCall(PetscInfo(tao,"Solving with %s formulation\n",TaoALMMTypes[auglag->type])); 37 PetscCall(TaoLogConvergenceHistory(tao, auglag->Lval, auglag->gnorm, auglag->cnorm, tao->ksp_its)); 38 PetscCall(TaoMonitor(tao, tao->niter, auglag->fval, auglag->gnorm, auglag->cnorm, 0.0)); 39 PetscCall((*tao->ops->convergencetest)(tao, tao->cnvP)); 40 /* set initial penalty factor and inner solver tolerance */ 41 switch (auglag->type) { 42 case TAO_ALMM_CLASSIC: 43 auglag->mu = auglag->mu0; 44 break; 45 case TAO_ALMM_PHR: 46 auglag->cenorm = 0.0; 47 if (tao->eq_constrained) { 48 PetscCall(VecDot(auglag->Ce, auglag->Ce, &auglag->cenorm)); 49 } 50 auglag->cinorm = 0.0; 51 if (tao->ineq_constrained) { 52 PetscCall(VecCopy(auglag->Ci, auglag->Ciwork)); 53 PetscCall(VecScale(auglag->Ciwork, -1.0)); 54 PetscCall(VecPointwiseMax(auglag->Ciwork, auglag->Cizero, auglag->Ciwork)); 55 PetscCall(VecDot(auglag->Ciwork, auglag->Ciwork, &auglag->cinorm)); 56 } 57 /* determine initial penalty factor based on the balance of constraint violation and objective function value */ 58 auglag->mu = PetscMax(1.e-6, PetscMin(10.0, 2.0*PetscAbsReal(auglag->fval)/(auglag->cenorm + auglag->cinorm))); 59 break; 60 default: 61 break; 62 } 63 auglag->gtol = auglag->gtol0; 64 PetscCall(PetscInfo(tao,"Initial penalty: %.2g\n",(double)auglag->mu)); 65 66 /* start aug-lag outer loop */ 67 while (tao->reason == TAO_CONTINUE_ITERATING) { 68 ++tao->niter; 69 /* update subsolver tolerance */ 70 PetscCall(PetscInfo(tao,"Subsolver tolerance: ||G|| <= %e\n",(double)auglag->gtol)); 71 PetscCall(TaoSetTolerances(auglag->subsolver, auglag->gtol, 0.0, 0.0)); 72 /* solve the bound-constrained or unconstrained subproblem */ 73 PetscCall(TaoSolve(auglag->subsolver)); 74 PetscCall(TaoGetConvergedReason(auglag->subsolver, &reason)); 75 tao->ksp_its += auglag->subsolver->ksp_its; 76 if (reason != TAO_CONVERGED_GATOL) { 77 PetscCall(PetscInfo(tao,"Subsolver failed to converge, reason: %s\n",TaoConvergedReasons[reason])); 78 } 79 /* evaluate solution and test convergence */ 80 PetscCall((*auglag->sub_obj)(tao)); 81 PetscCall(TaoALMMComputeOptimalityNorms_Private(tao)); 82 /* decide whether to update multipliers or not */ 83 updated = 0.0; 84 if (auglag->cnorm <= auglag->ytol) { 85 PetscCall(PetscInfo(tao,"Multipliers updated: ||C|| <= %e\n",(double)auglag->ytol)); 86 /* constraints are good, update multipliers and convergence tolerances */ 87 if (tao->eq_constrained) { 88 PetscCall(VecAXPY(auglag->Ye, auglag->mu, auglag->Ce)); 89 PetscCall(VecSet(auglag->Cework, auglag->ye_max)); 90 PetscCall(VecPointwiseMin(auglag->Ye, auglag->Cework, auglag->Ye)); 91 PetscCall(VecSet(auglag->Cework, auglag->ye_min)); 92 PetscCall(VecPointwiseMax(auglag->Ye, auglag->Cework, auglag->Ye)); 93 } 94 if (tao->ineq_constrained) { 95 PetscCall(VecAXPY(auglag->Yi, auglag->mu, auglag->Ci)); 96 PetscCall(VecSet(auglag->Ciwork, auglag->yi_max)); 97 PetscCall(VecPointwiseMin(auglag->Yi, auglag->Ciwork, auglag->Yi)); 98 PetscCall(VecSet(auglag->Ciwork, auglag->yi_min)); 99 PetscCall(VecPointwiseMax(auglag->Yi, auglag->Ciwork, auglag->Yi)); 100 } 101 /* tolerances are updated only for non-PHR methods */ 102 if (auglag->type != TAO_ALMM_PHR) { 103 auglag->ytol = PetscMax(tao->catol, auglag->ytol/PetscPowReal(auglag->mu, auglag->mu_pow_good)); 104 auglag->gtol = PetscMax(tao->gatol, auglag->gtol/auglag->mu); 105 } 106 updated = 1.0; 107 } else { 108 /* constraints are bad, update penalty factor */ 109 auglag->mu = PetscMin(auglag->mu_max, auglag->mu_fac*auglag->mu); 110 /* tolerances are reset only for non-PHR methods */ 111 if (auglag->type != TAO_ALMM_PHR) { 112 auglag->ytol = PetscMax(tao->catol, 0.1/PetscPowReal(auglag->mu, auglag->mu_pow_bad)); 113 auglag->gtol = PetscMax(tao->gatol, 1.0/auglag->mu); 114 } 115 PetscCall(PetscInfo(tao,"Penalty increased: mu = %.2g\n",(double)auglag->mu)); 116 } 117 PetscCall(TaoLogConvergenceHistory(tao, auglag->fval, auglag->gnorm, auglag->cnorm, tao->ksp_its)); 118 PetscCall(TaoMonitor(tao, tao->niter, auglag->fval, auglag->gnorm, auglag->cnorm, updated)); 119 PetscCall((*tao->ops->convergencetest)(tao, tao->cnvP)); 120 } 121 122 PetscFunctionReturn(0); 123 } 124 125 static PetscErrorCode TaoView_ALMM(Tao tao,PetscViewer viewer) 126 { 127 TAO_ALMM *auglag = (TAO_ALMM*)tao->data; 128 PetscBool isascii; 129 130 PetscFunctionBegin; 131 PetscCall(PetscViewerASCIIPushTab(viewer)); 132 PetscCall(TaoView(auglag->subsolver,viewer)); 133 PetscCall(PetscViewerASCIIPopTab(viewer)); 134 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 135 if (isascii) { 136 PetscCall(PetscViewerASCIIPushTab(viewer)); 137 PetscCall(PetscViewerASCIIPrintf(viewer, "ALMM Formulation Type: %s\n", TaoALMMTypes[auglag->type])); 138 PetscCall(PetscViewerASCIIPopTab(viewer)); 139 } 140 PetscFunctionReturn(0); 141 } 142 143 static PetscErrorCode TaoSetUp_ALMM(Tao tao) 144 { 145 TAO_ALMM *auglag = (TAO_ALMM*)tao->data; 146 VecType vec_type; 147 Vec SL, SU; 148 PetscBool is_cg = PETSC_FALSE, is_lmvm = PETSC_FALSE; 149 150 PetscFunctionBegin; 151 PetscCheck(!tao->ineq_doublesided,PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "TAOALMM does not support double-sided inequality constraint definition. Please restructure your inequality constrainst to fit the form c(x) >= 0."); 152 PetscCheck(tao->eq_constrained || tao->ineq_constrained,PetscObjectComm((PetscObject)tao), PETSC_ERR_ORDER, "Equality and/or inequality constraints must be defined before solver setup."); 153 PetscCall(TaoComputeVariableBounds(tao)); 154 /* alias base vectors and create extras */ 155 PetscCall(VecGetType(tao->solution, &vec_type)); 156 auglag->Px = tao->solution; 157 if (!tao->gradient) { /* base gradient */ 158 PetscCall(VecDuplicate(tao->solution, &tao->gradient)); 159 } 160 auglag->LgradX = tao->gradient; 161 if (!auglag->Xwork) { /* opt var work vector */ 162 PetscCall(VecDuplicate(tao->solution, &auglag->Xwork)); 163 } 164 if (tao->eq_constrained) { 165 auglag->Ce = tao->constraints_equality; 166 auglag->Ae = tao->jacobian_equality; 167 if (!auglag->Ye) { /* equality multipliers */ 168 PetscCall(VecDuplicate(auglag->Ce, &auglag->Ye)); 169 } 170 if (!auglag->Cework) { 171 PetscCall(VecDuplicate(auglag->Ce, &auglag->Cework)); 172 } 173 } 174 if (tao->ineq_constrained) { 175 auglag->Ci = tao->constraints_inequality; 176 auglag->Ai = tao->jacobian_inequality; 177 if (!auglag->Yi) { /* inequality multipliers */ 178 PetscCall(VecDuplicate(auglag->Ci, &auglag->Yi)); 179 } 180 if (!auglag->Ciwork) { 181 PetscCall(VecDuplicate(auglag->Ci, &auglag->Ciwork)); 182 } 183 if (!auglag->Cizero) { 184 PetscCall(VecDuplicate(auglag->Ci, &auglag->Cizero)); 185 PetscCall(VecZeroEntries(auglag->Cizero)); 186 } 187 if (!auglag->Ps) { /* slack vars */ 188 PetscCall(VecDuplicate(auglag->Ci, &auglag->Ps)); 189 } 190 if (!auglag->LgradS) { /* slack component of Lagrangian gradient */ 191 PetscCall(VecDuplicate(auglag->Ci, &auglag->LgradS)); 192 } 193 /* create vector for combined primal space and the associated communication objects */ 194 if (!auglag->P) { 195 PetscCall(PetscMalloc1(2, &auglag->Parr)); 196 auglag->Parr[0] = auglag->Px; auglag->Parr[1] = auglag->Ps; 197 PetscCall(VecConcatenate(2, auglag->Parr, &auglag->P, &auglag->Pis)); 198 PetscCall(PetscMalloc1(2, &auglag->Pscatter)); 199 PetscCall(VecScatterCreate(auglag->P, auglag->Pis[0], auglag->Px, NULL, &auglag->Pscatter[0])); 200 PetscCall(VecScatterCreate(auglag->P, auglag->Pis[1], auglag->Ps, NULL, &auglag->Pscatter[1])); 201 } 202 if (tao->eq_constrained) { 203 /* create vector for combined dual space and the associated communication objects */ 204 if (!auglag->Y) { 205 PetscCall(PetscMalloc1(2, &auglag->Yarr)); 206 auglag->Yarr[0] = auglag->Ye; auglag->Yarr[1] = auglag->Yi; 207 PetscCall(VecConcatenate(2, auglag->Yarr, &auglag->Y, &auglag->Yis)); 208 PetscCall(PetscMalloc1(2, &auglag->Yscatter)); 209 PetscCall(VecScatterCreate(auglag->Y, auglag->Yis[0], auglag->Ye, NULL, &auglag->Yscatter[0])); 210 PetscCall(VecScatterCreate(auglag->Y, auglag->Yis[1], auglag->Yi, NULL, &auglag->Yscatter[1])); 211 } 212 if (!auglag->C) { 213 PetscCall(VecDuplicate(auglag->Y, &auglag->C)); 214 } 215 } else { 216 if (!auglag->C) { 217 auglag->C = auglag->Ci; 218 } 219 if (!auglag->Y) { 220 auglag->Y = auglag->Yi; 221 } 222 } 223 } else { 224 if (!auglag->P) { 225 auglag->P = auglag->Px; 226 } 227 if (!auglag->G) { 228 auglag->G = auglag->LgradX; 229 } 230 if (!auglag->C) { 231 auglag->C = auglag->Ce; 232 } 233 if (!auglag->Y) { 234 auglag->Y = auglag->Ye; 235 } 236 } 237 /* initialize parameters */ 238 if (auglag->type == TAO_ALMM_PHR) { 239 auglag->mu_fac = 10.0; 240 auglag->yi_min = 0.0; 241 auglag->ytol0 = 0.5; 242 auglag->gtol0 = tao->gatol; 243 if (tao->gatol_changed && tao->catol_changed) { 244 PetscCall(PetscInfo(tao,"TAOALMM with PHR: different gradient and constraint tolerances are not supported, setting catol = gatol\n")); 245 tao->catol = tao->gatol; 246 } 247 } 248 /* set the Lagrangian formulation type for the subsolver */ 249 switch (auglag->type) { 250 case TAO_ALMM_CLASSIC: 251 auglag->sub_obj = TaoALMMComputeAugLagAndGradient_Private; 252 break; 253 case TAO_ALMM_PHR: 254 auglag->sub_obj = TaoALMMComputePHRLagAndGradient_Private; 255 break; 256 default: 257 break; 258 } 259 /* set up the subsolver */ 260 PetscCall(TaoSetSolution(auglag->subsolver, auglag->P)); 261 PetscCall(TaoSetObjective(auglag->subsolver, TaoALMMSubsolverObjective_Private, (void*)auglag)); 262 PetscCall(TaoSetObjectiveAndGradient(auglag->subsolver, NULL, TaoALMMSubsolverObjectiveAndGradient_Private, (void*)auglag)); 263 if (tao->bounded) { 264 /* make sure that the subsolver is a bound-constrained method */ 265 PetscCall(PetscObjectTypeCompare((PetscObject)auglag->subsolver, TAOCG, &is_cg)); 266 PetscCall(PetscObjectTypeCompare((PetscObject)auglag->subsolver, TAOLMVM, &is_lmvm)); 267 if (is_cg) { 268 PetscCall(TaoSetType(auglag->subsolver, TAOBNCG)); 269 PetscCall(PetscInfo(tao,"TAOCG detected for bound-constrained problem, switching to TAOBNCG instead.")); 270 } 271 if (is_lmvm) { 272 PetscCall(TaoSetType(auglag->subsolver, TAOBQNLS)); 273 PetscCall(PetscInfo(tao,"TAOLMVM detected for bound-constrained problem, switching to TAOBQNLS instead.")); 274 } 275 /* create lower and upper bound clone vectors for subsolver */ 276 if (!auglag->PL) { 277 PetscCall(VecDuplicate(auglag->P, &auglag->PL)); 278 } 279 if (!auglag->PU) { 280 PetscCall(VecDuplicate(auglag->P, &auglag->PU)); 281 } 282 if (tao->ineq_constrained) { 283 /* create lower and upper bounds for slack, set lower to 0 */ 284 PetscCall(VecDuplicate(auglag->Ci, &SL)); 285 PetscCall(VecSet(SL, 0.0)); 286 PetscCall(VecDuplicate(auglag->Ci, &SU)); 287 PetscCall(VecSet(SU, PETSC_INFINITY)); 288 /* combine opt var bounds with slack bounds */ 289 PetscCall(TaoALMMCombinePrimal_Private(tao, tao->XL, SL, auglag->PL)); 290 PetscCall(TaoALMMCombinePrimal_Private(tao, tao->XU, SU, auglag->PU)); 291 /* destroy work vectors */ 292 PetscCall(VecDestroy(&SL)); 293 PetscCall(VecDestroy(&SU)); 294 } else { 295 /* no inequality constraints, just copy bounds into the subsolver */ 296 PetscCall(VecCopy(tao->XL, auglag->PL)); 297 PetscCall(VecCopy(tao->XU, auglag->PU)); 298 } 299 PetscCall(TaoSetVariableBounds(auglag->subsolver, auglag->PL, auglag->PU)); 300 } 301 PetscCall(TaoSetUp(auglag->subsolver)); 302 auglag->G = auglag->subsolver->gradient; 303 304 PetscFunctionReturn(0); 305 } 306 307 static PetscErrorCode TaoDestroy_ALMM(Tao tao) 308 { 309 TAO_ALMM *auglag = (TAO_ALMM*)tao->data; 310 311 PetscFunctionBegin; 312 PetscCall(TaoDestroy(&auglag->subsolver)); 313 if (tao->setupcalled) { 314 PetscCall(VecDestroy(&auglag->Xwork)); /* opt work */ 315 if (tao->eq_constrained) { 316 PetscCall(VecDestroy(&auglag->Ye)); /* equality multipliers */ 317 PetscCall(VecDestroy(&auglag->Cework)); /* equality work vector */ 318 } 319 if (tao->ineq_constrained) { 320 PetscCall(VecDestroy(&auglag->Ps)); /* slack vars */ 321 auglag->Parr[0] = NULL; /* clear pointer to tao->solution, will be destroyed by TaoDestroy() shell */ 322 PetscCall(PetscFree(auglag->Parr)); /* array of primal vectors */ 323 PetscCall(VecDestroy(&auglag->LgradS)); /* slack grad */ 324 PetscCall(VecDestroy(&auglag->Cizero)); /* zero vector for pointwise max */ 325 PetscCall(VecDestroy(&auglag->Yi)); /* inequality multipliers */ 326 PetscCall(VecDestroy(&auglag->Ciwork)); /* inequality work vector */ 327 PetscCall(VecDestroy(&auglag->P)); /* combo primal */ 328 PetscCall(ISDestroy(&auglag->Pis[0])); /* index set for X inside P */ 329 PetscCall(ISDestroy(&auglag->Pis[1])); /* index set for S inside P */ 330 PetscCall(PetscFree(auglag->Pis)); /* array of P index sets */ 331 PetscCall(VecScatterDestroy(&auglag->Pscatter[0])); 332 PetscCall(VecScatterDestroy(&auglag->Pscatter[1])); 333 PetscCall(PetscFree(auglag->Pscatter)); 334 if (tao->eq_constrained) { 335 PetscCall(VecDestroy(&auglag->Y)); /* combo multipliers */ 336 PetscCall(PetscFree(auglag->Yarr)); /* array of dual vectors */ 337 PetscCall(VecDestroy(&auglag->C)); /* combo constraints */ 338 PetscCall(ISDestroy(&auglag->Yis[0])); /* index set for Ye inside Y */ 339 PetscCall(ISDestroy(&auglag->Yis[1])); /* index set for Yi inside Y */ 340 PetscCall(PetscFree(auglag->Yis)); 341 PetscCall(VecScatterDestroy(&auglag->Yscatter[0])); 342 PetscCall(VecScatterDestroy(&auglag->Yscatter[1])); 343 PetscCall(PetscFree(auglag->Yscatter)); 344 } 345 } 346 if (tao->bounded) { 347 PetscCall(VecDestroy(&auglag->PL)); /* lower bounds for subsolver */ 348 PetscCall(VecDestroy(&auglag->PU)); /* upper bounds for subsolver */ 349 } 350 } 351 PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoALMMGetType_C",NULL)); 352 PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoALMMSetType_C",NULL)); 353 PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoALMMGetSubsolver_C",NULL)); 354 PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoALMMSetSubsolver_C",NULL)); 355 PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoALMMGetMultipliers_C",NULL)); 356 PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoALMMSetMultipliers_C",NULL)); 357 PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoALMMGetPrimalIS_C",NULL)); 358 PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoALMMGetDualIS_C",NULL)); 359 PetscCall(PetscFree(tao->data)); 360 PetscFunctionReturn(0); 361 } 362 363 static PetscErrorCode TaoSetFromOptions_ALMM(PetscOptionItems *PetscOptionsObject,Tao tao) 364 { 365 TAO_ALMM *auglag = (TAO_ALMM*)tao->data; 366 PetscInt i; 367 368 PetscFunctionBegin; 369 PetscOptionsHeadBegin(PetscOptionsObject,"Augmented Lagrangian multiplier method solves problems with general constraints by converting them into a sequence of unconstrained problems."); 370 PetscCall(PetscOptionsReal("-tao_almm_mu_init","initial penalty parameter","",auglag->mu0,&auglag->mu0,NULL)); 371 PetscCall(PetscOptionsReal("-tao_almm_mu_factor","increase factor for the penalty parameter","",auglag->mu_fac,&auglag->mu_fac,NULL)); 372 PetscCall(PetscOptionsReal("-tao_almm_mu_power_good","exponential for penalty parameter when multiplier update is accepted","",auglag->mu_pow_good,&auglag->mu_pow_good,NULL)); 373 PetscCall(PetscOptionsReal("-tao_almm_mu_power_bad","exponential for penalty parameter when multiplier update is rejected","",auglag->mu_pow_bad,&auglag->mu_pow_bad,NULL)); 374 PetscCall(PetscOptionsReal("-tao_almm_mu_max","maximum safeguard for penalty parameter updates","",auglag->mu_max,&auglag->mu_max,NULL)); 375 PetscCall(PetscOptionsReal("-tao_almm_ye_min","minimum safeguard for equality multiplier updates","",auglag->ye_min,&auglag->ye_min,NULL)); 376 PetscCall(PetscOptionsReal("-tao_almm_ye_max","maximum safeguard for equality multipliers updates","",auglag->ye_max,&auglag->ye_max,NULL)); 377 PetscCall(PetscOptionsReal("-tao_almm_yi_min","minimum safeguard for inequality multipliers updates","",auglag->yi_min,&auglag->yi_min,NULL)); 378 PetscCall(PetscOptionsReal("-tao_almm_yi_max","maximum safeguard for inequality multipliers updates","",auglag->yi_max,&auglag->yi_max,NULL)); 379 PetscCall(PetscOptionsEnum("-tao_almm_type","augmented Lagrangian formulation type for the subproblem","TaoALMMType",TaoALMMTypes,(PetscEnum)auglag->type,(PetscEnum*)&auglag->type,NULL)); 380 PetscOptionsHeadEnd(); 381 PetscCall(TaoSetOptionsPrefix(auglag->subsolver,((PetscObject)tao)->prefix)); 382 PetscCall(TaoAppendOptionsPrefix(auglag->subsolver,"tao_almm_subsolver_")); 383 PetscCall(TaoSetFromOptions(auglag->subsolver)); 384 for (i=0; i<tao->numbermonitors; i++) { 385 PetscCall(PetscObjectReference((PetscObject)tao->monitorcontext[i])); 386 PetscCall(TaoSetMonitor(auglag->subsolver, tao->monitor[i], tao->monitorcontext[i], tao->monitordestroy[i])); 387 if (tao->monitor[i] == TaoMonitorDefault || tao->monitor[i] == TaoDefaultCMonitor || tao->monitor[i] == TaoDefaultGMonitor || tao->monitor[i] == TaoDefaultSMonitor) { 388 auglag->info = PETSC_TRUE; 389 } 390 } 391 PetscFunctionReturn(0); 392 } 393 394 /* -------------------------------------------------------- */ 395 396 /*MC 397 TaoALMM - Augmented Lagrangian multiplier method for solving nonlinear optimization problems with general constraints. 398 399 Options Database Keys: 400 + -tao_almm_mu_init <real> - initial penalty parameter (default: 10.) 401 . -tao_almm_mu_factor <real> - increase factor for the penalty parameter (default: 100.) 402 . -tao_almm_mu_max <real> - maximum safeguard for penalty parameter updates (default: 1.e20) 403 . -tao_almm_mu_power_good <real> - exponential for penalty parameter when multiplier update is accepted (default: 0.9) 404 . -tao_almm_mu_power_bad <real> - exponential for penalty parameter when multiplier update is rejected (default: 0.1) 405 . -tao_almm_ye_min <real> - minimum safeguard for equality multiplier updates (default: -1.e20) 406 . -tao_almm_ye_max <real> - maximum safeguard for equality multiplier updates (default: 1.e20) 407 . -tao_almm_yi_min <real> - minimum safeguard for inequality multiplier updates (default: -1.e20) 408 . -tao_almm_yi_max <real> - maximum safeguard for inequality multiplier updates (default: 1.e20) 409 - -tao_almm_type <classic,phr> - change formulation of the augmented Lagrangian merit function for the subproblem (default: classic) 410 411 Level: beginner 412 413 Notes: 414 This method converts a constrained problem into a sequence of unconstrained problems via the augmented 415 Lagrangian merit function. Bound constraints are pushed down to the subproblem without any modifications. 416 417 Two formulations are offered for the subproblem: canonical Hestenes-Powell augmented Lagrangian with slack 418 variables for inequality constraints, and a slack-less Powell-Hestenes-Rockafellar (PHR) formulation utilizing a 419 pointwise max() penalty on inequality constraints. The canonical augmented Lagrangian formulation typically 420 converges faster for most problems. However, PHR may be desirable for problems featuring a large number 421 of inequality constraints because it avoids inflating the size of the subproblem with slack variables. 422 423 The subproblem is solved using a nested first-order TAO solver. The user can retrieve a pointer to 424 the subsolver via TaoALMMGetSubsolver() or pass command line arguments to it using the 425 "-tao_almm_subsolver_" prefix. Currently, TaoALMM does not support second-order methods for the 426 subproblem. It is also highly recommended that the subsolver chosen by the user utilize a trust-region 427 strategy for globalization (default: TAOBQNKTR) especially if the outer problem features bound constraints. 428 429 .vb 430 while unconverged 431 solve argmin_x L(x) s.t. l <= x <= u 432 if ||c|| <= y_tol 433 if ||c|| <= c_tol && ||Lgrad|| <= g_tol: 434 problem converged, return solution 435 else 436 constraints sufficiently improved 437 update multipliers and tighten tolerances 438 endif 439 else 440 constraints did not improve 441 update penalty and loosen tolerances 442 endif 443 endwhile 444 .ve 445 446 .seealso: `TaoALMMGetType()`, `TaoALMMSetType()`, `TaoALMMSetSubsolver()`, `TaoALMMGetSubsolver()`, 447 `TaoALMMGetMultipliers()`, `TaoALMMSetMultipliers()`, `TaoALMMGetPrimalIS()`, `TaoALMMGetDualIS()` 448 M*/ 449 PETSC_EXTERN PetscErrorCode TaoCreate_ALMM(Tao tao) 450 { 451 TAO_ALMM *auglag; 452 453 PetscFunctionBegin; 454 PetscCall(PetscNewLog(tao, &auglag)); 455 456 tao->ops->destroy = TaoDestroy_ALMM; 457 tao->ops->setup = TaoSetUp_ALMM; 458 tao->ops->setfromoptions = TaoSetFromOptions_ALMM; 459 tao->ops->view = TaoView_ALMM; 460 tao->ops->solve = TaoSolve_ALMM; 461 462 tao->gatol = 1.e-5; 463 tao->grtol = 0.0; 464 tao->gttol = 0.0; 465 tao->catol = 1.e-5; 466 tao->crtol = 0.0; 467 468 tao->data = (void*)auglag; 469 auglag->parent = tao; 470 auglag->mu0 = 10.0; 471 auglag->mu = auglag->mu0; 472 auglag->mu_fac = 10.0; 473 auglag->mu_max = PETSC_INFINITY; 474 auglag->mu_pow_good = 0.9; 475 auglag->mu_pow_bad = 0.1; 476 auglag->ye_min = PETSC_NINFINITY; 477 auglag->ye_max = PETSC_INFINITY; 478 auglag->yi_min = PETSC_NINFINITY; 479 auglag->yi_max = PETSC_INFINITY; 480 auglag->ytol0 = 0.1/PetscPowReal(auglag->mu0, auglag->mu_pow_bad); 481 auglag->ytol = auglag->ytol0; 482 auglag->gtol0 = 1.0/auglag->mu0; 483 auglag->gtol = auglag->gtol0; 484 485 auglag->sub_obj = TaoALMMComputeAugLagAndGradient_Private; 486 auglag->type = TAO_ALMM_CLASSIC; 487 auglag->info = PETSC_FALSE; 488 489 PetscCall(TaoCreate(PetscObjectComm((PetscObject)tao),&auglag->subsolver)); 490 PetscCall(TaoSetType(auglag->subsolver, TAOBQNKTR)); 491 PetscCall(TaoSetTolerances(auglag->subsolver, auglag->gtol, 0.0, 0.0)); 492 PetscCall(TaoSetMaximumIterations(auglag->subsolver, 1000)); 493 PetscCall(TaoSetMaximumFunctionEvaluations(auglag->subsolver, 10000)); 494 PetscCall(TaoSetFunctionLowerBound(auglag->subsolver, PETSC_NINFINITY)); 495 PetscCall(PetscObjectIncrementTabLevel((PetscObject)auglag->subsolver,(PetscObject)tao,1)); 496 497 PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetType_C", TaoALMMGetType_Private)); 498 PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMSetType_C", TaoALMMSetType_Private)); 499 PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetSubsolver_C", TaoALMMGetSubsolver_Private)); 500 PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMSetSubsolver_C", TaoALMMSetSubsolver_Private)); 501 PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetMultipliers_C", TaoALMMGetMultipliers_Private)); 502 PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMSetMultipliers_C", TaoALMMSetMultipliers_Private)); 503 PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetPrimalIS_C", TaoALMMGetPrimalIS_Private)); 504 PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetDualIS_C", TaoALMMGetDualIS_Private)); 505 PetscFunctionReturn(0); 506 } 507 508 static PetscErrorCode TaoALMMCombinePrimal_Private(Tao tao, Vec X, Vec S, Vec P) 509 { 510 TAO_ALMM *auglag = (TAO_ALMM*)tao->data; 511 512 PetscFunctionBegin; 513 if (tao->ineq_constrained) { 514 PetscCall(VecScatterBegin(auglag->Pscatter[0], X, P, INSERT_VALUES, SCATTER_REVERSE)); 515 PetscCall(VecScatterEnd(auglag->Pscatter[0], X, P, INSERT_VALUES, SCATTER_REVERSE)); 516 PetscCall(VecScatterBegin(auglag->Pscatter[1], S, P, INSERT_VALUES, SCATTER_REVERSE)); 517 PetscCall(VecScatterEnd(auglag->Pscatter[1], S, P, INSERT_VALUES, SCATTER_REVERSE)); 518 } else { 519 PetscCall(VecCopy(X, P)); 520 } 521 PetscFunctionReturn(0); 522 } 523 524 static PetscErrorCode TaoALMMCombineDual_Private(Tao tao, Vec EQ, Vec IN, Vec Y) 525 { 526 TAO_ALMM *auglag = (TAO_ALMM*)tao->data; 527 528 PetscFunctionBegin; 529 if (tao->eq_constrained) { 530 if (tao->ineq_constrained) { 531 PetscCall(VecScatterBegin(auglag->Yscatter[0], EQ, Y, INSERT_VALUES, SCATTER_REVERSE)); 532 PetscCall(VecScatterEnd(auglag->Yscatter[0], EQ, Y, INSERT_VALUES, SCATTER_REVERSE)); 533 PetscCall(VecScatterBegin(auglag->Yscatter[1], IN, Y, INSERT_VALUES, SCATTER_REVERSE)); 534 PetscCall(VecScatterEnd(auglag->Yscatter[1], IN, Y, INSERT_VALUES, SCATTER_REVERSE)); 535 } else { 536 PetscCall(VecCopy(EQ, Y)); 537 } 538 } else { 539 PetscCall(VecCopy(IN, Y)); 540 } 541 PetscFunctionReturn(0); 542 } 543 544 static PetscErrorCode TaoALMMSplitPrimal_Private(Tao tao, Vec P, Vec X, Vec S) 545 { 546 TAO_ALMM *auglag = (TAO_ALMM*)tao->data; 547 548 PetscFunctionBegin; 549 if (tao->ineq_constrained) { 550 PetscCall(VecScatterBegin(auglag->Pscatter[0], P, X, INSERT_VALUES, SCATTER_FORWARD)); 551 PetscCall(VecScatterEnd(auglag->Pscatter[0], P, X, INSERT_VALUES, SCATTER_FORWARD)); 552 PetscCall(VecScatterBegin(auglag->Pscatter[1], P, S, INSERT_VALUES, SCATTER_FORWARD)); 553 PetscCall(VecScatterEnd(auglag->Pscatter[1], P, S, INSERT_VALUES, SCATTER_FORWARD)); 554 } else { 555 PetscCall(VecCopy(P, X)); 556 } 557 PetscFunctionReturn(0); 558 } 559 560 /* this assumes that the latest constraints are stored in Ce and Ci, and also combined in C */ 561 static PetscErrorCode TaoALMMComputeOptimalityNorms_Private(Tao tao) 562 { 563 TAO_ALMM *auglag = (TAO_ALMM*)tao->data; 564 565 PetscFunctionBegin; 566 /* if bounded, project the gradient */ 567 if (tao->bounded) { 568 PetscCall(VecBoundGradientProjection(auglag->LgradX, auglag->Px, tao->XL, tao->XU, auglag->LgradX)); 569 } 570 if (auglag->type == TAO_ALMM_PHR) { 571 PetscCall(VecNorm(auglag->LgradX, NORM_INFINITY, &auglag->gnorm)); 572 auglag->cenorm = 0.0; 573 if (tao->eq_constrained) { 574 PetscCall(VecNorm(auglag->Ce, NORM_INFINITY, &auglag->cenorm)); 575 } 576 auglag->cinorm = 0.0; 577 if (tao->ineq_constrained) { 578 PetscCall(VecCopy(auglag->Yi, auglag->Ciwork)); 579 PetscCall(VecScale(auglag->Ciwork, -1.0/auglag->mu)); 580 PetscCall(VecPointwiseMax(auglag->Ciwork, auglag->Ci, auglag->Ciwork)); 581 PetscCall(VecNorm(auglag->Ciwork, NORM_INFINITY, &auglag->cinorm)); 582 } 583 auglag->cnorm_old = auglag->cnorm; 584 auglag->cnorm = PetscMax(auglag->cenorm, auglag->cinorm); 585 auglag->ytol = auglag->ytol0 * auglag->cnorm_old; 586 } else { 587 PetscCall(VecNorm(auglag->LgradX, NORM_2, &auglag->gnorm)); 588 PetscCall(VecNorm(auglag->C, NORM_2, &auglag->cnorm)); 589 } 590 PetscFunctionReturn(0); 591 } 592 593 static PetscErrorCode TaoALMMEvaluateIterate_Private(Tao tao, Vec P) 594 { 595 TAO_ALMM *auglag = (TAO_ALMM*)tao->data; 596 597 PetscFunctionBegin; 598 /* split solution into primal and slack components */ 599 PetscCall(TaoALMMSplitPrimal_Private(tao, auglag->P, auglag->Px, auglag->Ps)); 600 601 /* compute f, df/dx and the constraints */ 602 PetscCall(TaoComputeObjectiveAndGradient(tao, auglag->Px, &auglag->fval, auglag->LgradX)); 603 if (tao->eq_constrained) { 604 PetscCall(TaoComputeEqualityConstraints(tao, auglag->Px, auglag->Ce)); 605 PetscCall(TaoComputeJacobianEquality(tao, auglag->Px, auglag->Ae, auglag->Ae)); 606 } 607 if (tao->ineq_constrained) { 608 PetscCall(TaoComputeInequalityConstraints(tao, auglag->Px, auglag->Ci)); 609 PetscCall(TaoComputeJacobianInequality(tao, auglag->Px, auglag->Ai, auglag->Ai)); 610 switch (auglag->type) { 611 case TAO_ALMM_CLASSIC: 612 /* classic formulation converts inequality to equality constraints via slack variables */ 613 PetscCall(VecAXPY(auglag->Ci, -1.0, auglag->Ps)); 614 break; 615 case TAO_ALMM_PHR: 616 /* PHR is based on Ci <= 0 while TAO defines Ci >= 0 so we hit it with a negative sign */ 617 PetscCall(VecScale(auglag->Ci, -1.0)); 618 PetscCall(MatScale(auglag->Ai, -1.0)); 619 break; 620 default: 621 break; 622 } 623 } 624 /* combine constraints into one vector */ 625 PetscCall(TaoALMMCombineDual_Private(tao, auglag->Ce, auglag->Ci, auglag->C)); 626 PetscFunctionReturn(0); 627 } 628 629 /* 630 Lphr = f + 0.5*mu*[ (Ce + Ye/mu)^T (Ce + Ye/mu) + pmin(0, Ci + Yi/mu)^T pmin(0, Ci + Yi/mu)] 631 632 dLphr/dX = dF/dX + mu*[ (Ce + Ye/mu)^T Ae + pmin(0, Ci + Yi/mu)^T Ai] 633 634 dLphr/dS = 0 635 */ 636 static PetscErrorCode TaoALMMComputePHRLagAndGradient_Private(Tao tao) 637 { 638 TAO_ALMM *auglag = (TAO_ALMM*)tao->data; 639 PetscReal eq_norm=0.0, ineq_norm=0.0; 640 641 PetscFunctionBegin; 642 PetscCall(TaoALMMEvaluateIterate_Private(tao, auglag->P)); 643 if (tao->eq_constrained) { 644 /* Ce_work = mu*(Ce + Ye/mu) */ 645 PetscCall(VecWAXPY(auglag->Cework, 1.0/auglag->mu, auglag->Ye, auglag->Ce)); 646 PetscCall(VecDot(auglag->Cework, auglag->Cework, &eq_norm)); /* contribution to scalar Lagrangian */ 647 PetscCall(VecScale(auglag->Cework, auglag->mu)); 648 /* dL/dX += mu*(Ce + Ye/mu)^T Ae */ 649 PetscCall(MatMultTransposeAdd(auglag->Ae, auglag->Cework, auglag->LgradX, auglag->LgradX)); 650 } 651 if (tao->ineq_constrained) { 652 /* Ci_work = mu * pmax(0, Ci + Yi/mu) where pmax() is pointwise max() */ 653 PetscCall(VecWAXPY(auglag->Ciwork, 1.0/auglag->mu, auglag->Yi, auglag->Ci)); 654 PetscCall(VecPointwiseMax(auglag->Ciwork, auglag->Cizero, auglag->Ciwork)); 655 PetscCall(VecDot(auglag->Ciwork, auglag->Ciwork, &ineq_norm)); /* contribution to scalar Lagrangian */ 656 /* dL/dX += mu * pmax(0, Ci + Yi/mu)^T Ai */ 657 PetscCall(VecScale(auglag->Ciwork, auglag->mu)); 658 PetscCall(MatMultTransposeAdd(auglag->Ai, auglag->Ciwork, auglag->LgradX, auglag->LgradX)); 659 /* dL/dS = 0 because there are no slacks in PHR */ 660 PetscCall(VecZeroEntries(auglag->LgradS)); 661 } 662 /* combine gradient together */ 663 PetscCall(TaoALMMCombinePrimal_Private(tao, auglag->LgradX, auglag->LgradS, auglag->G)); 664 /* compute L = f + 0.5 * mu * [(Ce + Ye/mu)^T (Ce + Ye/mu) + pmax(0, Ci + Yi/mu)^T pmax(0, Ci + Yi/mu)] */ 665 auglag->Lval = auglag->fval + 0.5*auglag->mu*(eq_norm + ineq_norm); 666 PetscFunctionReturn(0); 667 } 668 669 /* 670 Lc = F + Ye^TCe + Yi^T(Ci - S) + 0.5*mu*[Ce^TCe + (Ci - S)^T(Ci - S)] 671 672 dLc/dX = dF/dX + Ye^TAe + Yi^TAi + 0.5*mu*[Ce^TAe + (Ci - S)^TAi] 673 674 dLc/dS = -[Yi + mu*(Ci - S)] 675 */ 676 static PetscErrorCode TaoALMMComputeAugLagAndGradient_Private(Tao tao) 677 { 678 TAO_ALMM *auglag = (TAO_ALMM*)tao->data; 679 PetscReal yeTce=0.0, yiTcims=0.0, ceTce=0.0, cimsTcims=0.0; 680 681 PetscFunctionBegin; 682 PetscCall(TaoALMMEvaluateIterate_Private(tao, auglag->P)); 683 if (tao->eq_constrained) { 684 /* compute scalar contributions */ 685 PetscCall(VecDot(auglag->Ye, auglag->Ce, &yeTce)); 686 PetscCall(VecDot(auglag->Ce, auglag->Ce, &ceTce)); 687 /* dL/dX += ye^T Ae */ 688 PetscCall(MatMultTransposeAdd(auglag->Ae, auglag->Ye, auglag->LgradX, auglag->LgradX)); 689 /* dL/dX += mu * ce^T Ae */ 690 PetscCall(MatMultTranspose(auglag->Ae, auglag->Ce, auglag->Xwork)); 691 PetscCall(VecAXPY(auglag->LgradX, auglag->mu, auglag->Xwork)); 692 } 693 if (tao->ineq_constrained) { 694 /* compute scalar contributions */ 695 PetscCall(VecDot(auglag->Yi, auglag->Ci, &yiTcims)); 696 PetscCall(VecDot(auglag->Ci, auglag->Ci, &cimsTcims)); 697 /* dL/dX += yi^T Ai */ 698 PetscCall(MatMultTransposeAdd(auglag->Ai, auglag->Yi, auglag->LgradX, auglag->LgradX)); 699 /* dL/dX += mu * (ci - s)^T Ai */ 700 PetscCall(MatMultTranspose(auglag->Ai, auglag->Ci, auglag->Xwork)); 701 PetscCall(VecAXPY(auglag->LgradX, auglag->mu, auglag->Xwork)); 702 /* dL/dS = -[yi + mu*(ci - s)] */ 703 PetscCall(VecWAXPY(auglag->LgradS, auglag->mu, auglag->Ci, auglag->Yi)); 704 PetscCall(VecScale(auglag->LgradS, -1.0)); 705 } 706 /* combine gradient together */ 707 PetscCall(TaoALMMCombinePrimal_Private(tao, auglag->LgradX, auglag->LgradS, auglag->G)); 708 /* compute L = f + ye^T ce + yi^T (ci - s) + 0.5*mu*||ce||^2 + 0.5*mu*||ci - s||^2 */ 709 auglag->Lval = auglag->fval + yeTce + yiTcims + 0.5*auglag->mu*(ceTce + cimsTcims); 710 PetscFunctionReturn(0); 711 } 712 713 PetscErrorCode TaoALMMSubsolverObjective_Private(Tao tao, Vec P, PetscReal *Lval, void *ctx) 714 { 715 TAO_ALMM *auglag = (TAO_ALMM*)ctx; 716 717 PetscFunctionBegin; 718 PetscCall(VecCopy(P, auglag->P)); 719 PetscCall((*auglag->sub_obj)(auglag->parent)); 720 *Lval = auglag->Lval; 721 PetscFunctionReturn(0); 722 } 723 724 PetscErrorCode TaoALMMSubsolverObjectiveAndGradient_Private(Tao tao, Vec P, PetscReal *Lval, Vec G, void *ctx) 725 { 726 TAO_ALMM *auglag = (TAO_ALMM*)ctx; 727 728 PetscFunctionBegin; 729 PetscCall(VecCopy(P, auglag->P)); 730 PetscCall((*auglag->sub_obj)(auglag->parent)); 731 PetscCall(VecCopy(auglag->G, G)); 732 *Lval = auglag->Lval; 733 PetscFunctionReturn(0); 734 } 735