1661095bbSAlp Dener #include <../src/tao/constrained/impls/almm/almm.h> /*I "petsctao.h" I*/ /*I "petscvec.h" I*/ 2661095bbSAlp Dener #include <petsctao.h> 3661095bbSAlp Dener #include <petsc/private/petscimpl.h> 4661095bbSAlp Dener #include <petsc/private/vecimpl.h> 5661095bbSAlp Dener 6661095bbSAlp Dener static PetscErrorCode TaoALMMCombinePrimal_Private(Tao, Vec, Vec, Vec); 7661095bbSAlp Dener static PetscErrorCode TaoALMMCombineDual_Private(Tao, Vec, Vec, Vec); 8661095bbSAlp Dener static PetscErrorCode TaoALMMSplitPrimal_Private(Tao, Vec, Vec, Vec); 9661095bbSAlp Dener static PetscErrorCode TaoALMMComputeOptimalityNorms_Private(Tao); 10661095bbSAlp Dener static PetscErrorCode TaoALMMComputeAugLagAndGradient_Private(Tao); 11661095bbSAlp Dener static PetscErrorCode TaoALMMComputePHRLagAndGradient_Private(Tao); 12661095bbSAlp Dener 13d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSolve_ALMM(Tao tao) 14d71ae5a4SJacob Faibussowitsch { 15661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM *)tao->data; 16661095bbSAlp Dener TaoConvergedReason reason; 17661095bbSAlp Dener PetscReal updated; 18661095bbSAlp Dener 19661095bbSAlp Dener PetscFunctionBegin; 20661095bbSAlp Dener /* reset initial multiplier/slack guess */ 21661095bbSAlp Dener if (!tao->recycle) { 22661095bbSAlp Dener if (tao->ineq_constrained) { 239566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(auglag->Ps)); 249566063dSJacob Faibussowitsch PetscCall(TaoALMMCombinePrimal_Private(tao, auglag->Px, auglag->Ps, auglag->P)); 259566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(auglag->Yi)); 26661095bbSAlp Dener } 271baa6e33SBarry Smith if (tao->eq_constrained) PetscCall(VecZeroEntries(auglag->Ye)); 28661095bbSAlp Dener } 29661095bbSAlp Dener 30661095bbSAlp Dener /* compute initial nonlinear Lagrangian and its derivatives */ 319566063dSJacob Faibussowitsch PetscCall((*auglag->sub_obj)(tao)); 329566063dSJacob Faibussowitsch PetscCall(TaoALMMComputeOptimalityNorms_Private(tao)); 33661095bbSAlp Dener /* print initial step and check convergence */ 349566063dSJacob Faibussowitsch PetscCall(PetscInfo(tao, "Solving with %s formulation\n", TaoALMMTypes[auglag->type])); 359566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao, auglag->Lval, auglag->gnorm, auglag->cnorm, tao->ksp_its)); 369566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao, tao->niter, auglag->fval, auglag->gnorm, auglag->cnorm, 0.0)); 37dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, convergencetest, tao->cnvP); 38661095bbSAlp Dener /* set initial penalty factor and inner solver tolerance */ 39661095bbSAlp Dener switch (auglag->type) { 40d71ae5a4SJacob Faibussowitsch case TAO_ALMM_CLASSIC: 41d71ae5a4SJacob Faibussowitsch auglag->mu = auglag->mu0; 42d71ae5a4SJacob Faibussowitsch break; 437721a36fSStefano Zampini case TAO_ALMM_PHR: 44661095bbSAlp Dener auglag->cenorm = 0.0; 4548a46eb9SPierre Jolivet if (tao->eq_constrained) PetscCall(VecDot(auglag->Ce, auglag->Ce, &auglag->cenorm)); 46661095bbSAlp Dener auglag->cinorm = 0.0; 47661095bbSAlp Dener if (tao->ineq_constrained) { 489566063dSJacob Faibussowitsch PetscCall(VecCopy(auglag->Ci, auglag->Ciwork)); 499566063dSJacob Faibussowitsch PetscCall(VecScale(auglag->Ciwork, -1.0)); 509566063dSJacob Faibussowitsch PetscCall(VecPointwiseMax(auglag->Ciwork, auglag->Cizero, auglag->Ciwork)); 519566063dSJacob Faibussowitsch PetscCall(VecDot(auglag->Ciwork, auglag->Ciwork, &auglag->cinorm)); 52661095bbSAlp Dener } 53661095bbSAlp Dener /* determine initial penalty factor based on the balance of constraint violation and objective function value */ 54661095bbSAlp Dener auglag->mu = PetscMax(1.e-6, PetscMin(10.0, 2.0 * PetscAbsReal(auglag->fval) / (auglag->cenorm + auglag->cinorm))); 55661095bbSAlp Dener break; 56d71ae5a4SJacob Faibussowitsch default: 57d71ae5a4SJacob Faibussowitsch break; 58661095bbSAlp Dener } 59661095bbSAlp Dener auglag->gtol = auglag->gtol0; 6063a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(tao, "Initial penalty: %.2g\n", (double)auglag->mu)); 61661095bbSAlp Dener 62661095bbSAlp Dener /* start aug-lag outer loop */ 63661095bbSAlp Dener while (tao->reason == TAO_CONTINUE_ITERATING) { 64661095bbSAlp Dener ++tao->niter; 65661095bbSAlp Dener /* update subsolver tolerance */ 6663a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(tao, "Subsolver tolerance: ||G|| <= %e\n", (double)auglag->gtol)); 679566063dSJacob Faibussowitsch PetscCall(TaoSetTolerances(auglag->subsolver, auglag->gtol, 0.0, 0.0)); 68661095bbSAlp Dener /* solve the bound-constrained or unconstrained subproblem */ 699566063dSJacob Faibussowitsch PetscCall(TaoSolve(auglag->subsolver)); 709566063dSJacob Faibussowitsch PetscCall(TaoGetConvergedReason(auglag->subsolver, &reason)); 71661095bbSAlp Dener tao->ksp_its += auglag->subsolver->ksp_its; 7248a46eb9SPierre Jolivet if (reason != TAO_CONVERGED_GATOL) PetscCall(PetscInfo(tao, "Subsolver failed to converge, reason: %s\n", TaoConvergedReasons[reason])); 73661095bbSAlp Dener /* evaluate solution and test convergence */ 749566063dSJacob Faibussowitsch PetscCall((*auglag->sub_obj)(tao)); 759566063dSJacob Faibussowitsch PetscCall(TaoALMMComputeOptimalityNorms_Private(tao)); 76661095bbSAlp Dener /* decide whether to update multipliers or not */ 77661095bbSAlp Dener updated = 0.0; 78661095bbSAlp Dener if (auglag->cnorm <= auglag->ytol) { 7963a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(tao, "Multipliers updated: ||C|| <= %e\n", (double)auglag->ytol)); 80661095bbSAlp Dener /* constraints are good, update multipliers and convergence tolerances */ 81661095bbSAlp Dener if (tao->eq_constrained) { 829566063dSJacob Faibussowitsch PetscCall(VecAXPY(auglag->Ye, auglag->mu, auglag->Ce)); 839566063dSJacob Faibussowitsch PetscCall(VecSet(auglag->Cework, auglag->ye_max)); 849566063dSJacob Faibussowitsch PetscCall(VecPointwiseMin(auglag->Ye, auglag->Cework, auglag->Ye)); 859566063dSJacob Faibussowitsch PetscCall(VecSet(auglag->Cework, auglag->ye_min)); 869566063dSJacob Faibussowitsch PetscCall(VecPointwiseMax(auglag->Ye, auglag->Cework, auglag->Ye)); 87661095bbSAlp Dener } 88661095bbSAlp Dener if (tao->ineq_constrained) { 899566063dSJacob Faibussowitsch PetscCall(VecAXPY(auglag->Yi, auglag->mu, auglag->Ci)); 909566063dSJacob Faibussowitsch PetscCall(VecSet(auglag->Ciwork, auglag->yi_max)); 919566063dSJacob Faibussowitsch PetscCall(VecPointwiseMin(auglag->Yi, auglag->Ciwork, auglag->Yi)); 929566063dSJacob Faibussowitsch PetscCall(VecSet(auglag->Ciwork, auglag->yi_min)); 939566063dSJacob Faibussowitsch PetscCall(VecPointwiseMax(auglag->Yi, auglag->Ciwork, auglag->Yi)); 94661095bbSAlp Dener } 95661095bbSAlp Dener /* tolerances are updated only for non-PHR methods */ 96661095bbSAlp Dener if (auglag->type != TAO_ALMM_PHR) { 97661095bbSAlp Dener auglag->ytol = PetscMax(tao->catol, auglag->ytol / PetscPowReal(auglag->mu, auglag->mu_pow_good)); 98661095bbSAlp Dener auglag->gtol = PetscMax(tao->gatol, auglag->gtol / auglag->mu); 99661095bbSAlp Dener } 100661095bbSAlp Dener updated = 1.0; 101661095bbSAlp Dener } else { 102661095bbSAlp Dener /* constraints are bad, update penalty factor */ 103661095bbSAlp Dener auglag->mu = PetscMin(auglag->mu_max, auglag->mu_fac * auglag->mu); 104661095bbSAlp Dener /* tolerances are reset only for non-PHR methods */ 105661095bbSAlp Dener if (auglag->type != TAO_ALMM_PHR) { 106661095bbSAlp Dener auglag->ytol = PetscMax(tao->catol, 0.1 / PetscPowReal(auglag->mu, auglag->mu_pow_bad)); 107661095bbSAlp Dener auglag->gtol = PetscMax(tao->gatol, 1.0 / auglag->mu); 108661095bbSAlp Dener } 10963a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(tao, "Penalty increased: mu = %.2g\n", (double)auglag->mu)); 110661095bbSAlp Dener } 1119566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao, auglag->fval, auglag->gnorm, auglag->cnorm, tao->ksp_its)); 1129566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao, tao->niter, auglag->fval, auglag->gnorm, auglag->cnorm, updated)); 113dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, convergencetest, tao->cnvP); 114661095bbSAlp Dener } 115661095bbSAlp Dener 116661095bbSAlp Dener PetscFunctionReturn(0); 117661095bbSAlp Dener } 118661095bbSAlp Dener 119d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoView_ALMM(Tao tao, PetscViewer viewer) 120d71ae5a4SJacob Faibussowitsch { 121661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM *)tao->data; 122661095bbSAlp Dener PetscBool isascii; 123661095bbSAlp Dener 124661095bbSAlp Dener PetscFunctionBegin; 1259566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 1269566063dSJacob Faibussowitsch PetscCall(TaoView(auglag->subsolver, viewer)); 1279566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 1289566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 129661095bbSAlp Dener if (isascii) { 1309566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 1319566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "ALMM Formulation Type: %s\n", TaoALMMTypes[auglag->type])); 1329566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 133661095bbSAlp Dener } 134661095bbSAlp Dener PetscFunctionReturn(0); 135661095bbSAlp Dener } 136661095bbSAlp Dener 137d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSetUp_ALMM(Tao tao) 138d71ae5a4SJacob Faibussowitsch { 139661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM *)tao->data; 140661095bbSAlp Dener VecType vec_type; 141661095bbSAlp Dener Vec SL, SU; 142661095bbSAlp Dener PetscBool is_cg = PETSC_FALSE, is_lmvm = PETSC_FALSE; 143661095bbSAlp Dener 144661095bbSAlp Dener PetscFunctionBegin; 1453c859ba3SBarry Smith 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."); 1463c859ba3SBarry Smith PetscCheck(tao->eq_constrained || tao->ineq_constrained, PetscObjectComm((PetscObject)tao), PETSC_ERR_ORDER, "Equality and/or inequality constraints must be defined before solver setup."); 1479566063dSJacob Faibussowitsch PetscCall(TaoComputeVariableBounds(tao)); 148661095bbSAlp Dener /* alias base vectors and create extras */ 1499566063dSJacob Faibussowitsch PetscCall(VecGetType(tao->solution, &vec_type)); 150661095bbSAlp Dener auglag->Px = tao->solution; 151661095bbSAlp Dener if (!tao->gradient) { /* base gradient */ 1529566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &tao->gradient)); 153661095bbSAlp Dener } 154661095bbSAlp Dener auglag->LgradX = tao->gradient; 155661095bbSAlp Dener if (!auglag->Xwork) { /* opt var work vector */ 1569566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &auglag->Xwork)); 157661095bbSAlp Dener } 158661095bbSAlp Dener if (tao->eq_constrained) { 159661095bbSAlp Dener auglag->Ce = tao->constraints_equality; 160661095bbSAlp Dener auglag->Ae = tao->jacobian_equality; 161661095bbSAlp Dener if (!auglag->Ye) { /* equality multipliers */ 1629566063dSJacob Faibussowitsch PetscCall(VecDuplicate(auglag->Ce, &auglag->Ye)); 163661095bbSAlp Dener } 16448a46eb9SPierre Jolivet if (!auglag->Cework) PetscCall(VecDuplicate(auglag->Ce, &auglag->Cework)); 165661095bbSAlp Dener } 166661095bbSAlp Dener if (tao->ineq_constrained) { 167661095bbSAlp Dener auglag->Ci = tao->constraints_inequality; 168661095bbSAlp Dener auglag->Ai = tao->jacobian_inequality; 169661095bbSAlp Dener if (!auglag->Yi) { /* inequality multipliers */ 1709566063dSJacob Faibussowitsch PetscCall(VecDuplicate(auglag->Ci, &auglag->Yi)); 171661095bbSAlp Dener } 17248a46eb9SPierre Jolivet if (!auglag->Ciwork) PetscCall(VecDuplicate(auglag->Ci, &auglag->Ciwork)); 173661095bbSAlp Dener if (!auglag->Cizero) { 1749566063dSJacob Faibussowitsch PetscCall(VecDuplicate(auglag->Ci, &auglag->Cizero)); 1759566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(auglag->Cizero)); 176661095bbSAlp Dener } 177661095bbSAlp Dener if (!auglag->Ps) { /* slack vars */ 1789566063dSJacob Faibussowitsch PetscCall(VecDuplicate(auglag->Ci, &auglag->Ps)); 179661095bbSAlp Dener } 180661095bbSAlp Dener if (!auglag->LgradS) { /* slack component of Lagrangian gradient */ 1819566063dSJacob Faibussowitsch PetscCall(VecDuplicate(auglag->Ci, &auglag->LgradS)); 182661095bbSAlp Dener } 183661095bbSAlp Dener /* create vector for combined primal space and the associated communication objects */ 184661095bbSAlp Dener if (!auglag->P) { 1859566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2, &auglag->Parr)); 1869371c9d4SSatish Balay auglag->Parr[0] = auglag->Px; 1879371c9d4SSatish Balay auglag->Parr[1] = auglag->Ps; 1889566063dSJacob Faibussowitsch PetscCall(VecConcatenate(2, auglag->Parr, &auglag->P, &auglag->Pis)); 1899566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2, &auglag->Pscatter)); 1909566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(auglag->P, auglag->Pis[0], auglag->Px, NULL, &auglag->Pscatter[0])); 1919566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(auglag->P, auglag->Pis[1], auglag->Ps, NULL, &auglag->Pscatter[1])); 192661095bbSAlp Dener } 193661095bbSAlp Dener if (tao->eq_constrained) { 194661095bbSAlp Dener /* create vector for combined dual space and the associated communication objects */ 195661095bbSAlp Dener if (!auglag->Y) { 1969566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2, &auglag->Yarr)); 1979371c9d4SSatish Balay auglag->Yarr[0] = auglag->Ye; 1989371c9d4SSatish Balay auglag->Yarr[1] = auglag->Yi; 1999566063dSJacob Faibussowitsch PetscCall(VecConcatenate(2, auglag->Yarr, &auglag->Y, &auglag->Yis)); 2009566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2, &auglag->Yscatter)); 2019566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(auglag->Y, auglag->Yis[0], auglag->Ye, NULL, &auglag->Yscatter[0])); 2029566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(auglag->Y, auglag->Yis[1], auglag->Yi, NULL, &auglag->Yscatter[1])); 203661095bbSAlp Dener } 20448a46eb9SPierre Jolivet if (!auglag->C) PetscCall(VecDuplicate(auglag->Y, &auglag->C)); 2059371c9d4SSatish Balay } else { 206ad540459SPierre Jolivet if (!auglag->C) auglag->C = auglag->Ci; 207ad540459SPierre Jolivet if (!auglag->Y) auglag->Y = auglag->Yi; 208661095bbSAlp Dener } 209661095bbSAlp Dener } else { 210ad540459SPierre Jolivet if (!auglag->P) auglag->P = auglag->Px; 211ad540459SPierre Jolivet if (!auglag->G) auglag->G = auglag->LgradX; 212ad540459SPierre Jolivet if (!auglag->C) auglag->C = auglag->Ce; 213ad540459SPierre Jolivet if (!auglag->Y) auglag->Y = auglag->Ye; 214661095bbSAlp Dener } 215661095bbSAlp Dener /* initialize parameters */ 216661095bbSAlp Dener if (auglag->type == TAO_ALMM_PHR) { 217661095bbSAlp Dener auglag->mu_fac = 10.0; 218661095bbSAlp Dener auglag->yi_min = 0.0; 219661095bbSAlp Dener auglag->ytol0 = 0.5; 220661095bbSAlp Dener auglag->gtol0 = tao->gatol; 221661095bbSAlp Dener if (tao->gatol_changed && tao->catol_changed) { 2229566063dSJacob Faibussowitsch PetscCall(PetscInfo(tao, "TAOALMM with PHR: different gradient and constraint tolerances are not supported, setting catol = gatol\n")); 223661095bbSAlp Dener tao->catol = tao->gatol; 224661095bbSAlp Dener } 225661095bbSAlp Dener } 226661095bbSAlp Dener /* set the Lagrangian formulation type for the subsolver */ 227661095bbSAlp Dener switch (auglag->type) { 228d71ae5a4SJacob Faibussowitsch case TAO_ALMM_CLASSIC: 229d71ae5a4SJacob Faibussowitsch auglag->sub_obj = TaoALMMComputeAugLagAndGradient_Private; 230d71ae5a4SJacob Faibussowitsch break; 231d71ae5a4SJacob Faibussowitsch case TAO_ALMM_PHR: 232d71ae5a4SJacob Faibussowitsch auglag->sub_obj = TaoALMMComputePHRLagAndGradient_Private; 233d71ae5a4SJacob Faibussowitsch break; 234d71ae5a4SJacob Faibussowitsch default: 235d71ae5a4SJacob Faibussowitsch break; 236661095bbSAlp Dener } 237661095bbSAlp Dener /* set up the subsolver */ 2389566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(auglag->subsolver, auglag->P)); 2399566063dSJacob Faibussowitsch PetscCall(TaoSetObjective(auglag->subsolver, TaoALMMSubsolverObjective_Private, (void *)auglag)); 2409566063dSJacob Faibussowitsch PetscCall(TaoSetObjectiveAndGradient(auglag->subsolver, NULL, TaoALMMSubsolverObjectiveAndGradient_Private, (void *)auglag)); 241661095bbSAlp Dener if (tao->bounded) { 242661095bbSAlp Dener /* make sure that the subsolver is a bound-constrained method */ 2439566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)auglag->subsolver, TAOCG, &is_cg)); 2449566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)auglag->subsolver, TAOLMVM, &is_lmvm)); 245661095bbSAlp Dener if (is_cg) { 2469566063dSJacob Faibussowitsch PetscCall(TaoSetType(auglag->subsolver, TAOBNCG)); 2479566063dSJacob Faibussowitsch PetscCall(PetscInfo(tao, "TAOCG detected for bound-constrained problem, switching to TAOBNCG instead.")); 248661095bbSAlp Dener } 249661095bbSAlp Dener if (is_lmvm) { 2509566063dSJacob Faibussowitsch PetscCall(TaoSetType(auglag->subsolver, TAOBQNLS)); 2519566063dSJacob Faibussowitsch PetscCall(PetscInfo(tao, "TAOLMVM detected for bound-constrained problem, switching to TAOBQNLS instead.")); 252661095bbSAlp Dener } 253661095bbSAlp Dener /* create lower and upper bound clone vectors for subsolver */ 25448a46eb9SPierre Jolivet if (!auglag->PL) PetscCall(VecDuplicate(auglag->P, &auglag->PL)); 25548a46eb9SPierre Jolivet if (!auglag->PU) PetscCall(VecDuplicate(auglag->P, &auglag->PU)); 256661095bbSAlp Dener if (tao->ineq_constrained) { 257661095bbSAlp Dener /* create lower and upper bounds for slack, set lower to 0 */ 2589566063dSJacob Faibussowitsch PetscCall(VecDuplicate(auglag->Ci, &SL)); 2599566063dSJacob Faibussowitsch PetscCall(VecSet(SL, 0.0)); 2609566063dSJacob Faibussowitsch PetscCall(VecDuplicate(auglag->Ci, &SU)); 2619566063dSJacob Faibussowitsch PetscCall(VecSet(SU, PETSC_INFINITY)); 262661095bbSAlp Dener /* combine opt var bounds with slack bounds */ 2639566063dSJacob Faibussowitsch PetscCall(TaoALMMCombinePrimal_Private(tao, tao->XL, SL, auglag->PL)); 2649566063dSJacob Faibussowitsch PetscCall(TaoALMMCombinePrimal_Private(tao, tao->XU, SU, auglag->PU)); 265661095bbSAlp Dener /* destroy work vectors */ 2669566063dSJacob Faibussowitsch PetscCall(VecDestroy(&SL)); 2679566063dSJacob Faibussowitsch PetscCall(VecDestroy(&SU)); 268661095bbSAlp Dener } else { 269661095bbSAlp Dener /* no inequality constraints, just copy bounds into the subsolver */ 2709566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->XL, auglag->PL)); 2719566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->XU, auglag->PU)); 272661095bbSAlp Dener } 2739566063dSJacob Faibussowitsch PetscCall(TaoSetVariableBounds(auglag->subsolver, auglag->PL, auglag->PU)); 274661095bbSAlp Dener } 2759566063dSJacob Faibussowitsch PetscCall(TaoSetUp(auglag->subsolver)); 276661095bbSAlp Dener auglag->G = auglag->subsolver->gradient; 277661095bbSAlp Dener 278661095bbSAlp Dener PetscFunctionReturn(0); 279661095bbSAlp Dener } 280661095bbSAlp Dener 281d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoDestroy_ALMM(Tao tao) 282d71ae5a4SJacob Faibussowitsch { 283661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM *)tao->data; 284661095bbSAlp Dener 285661095bbSAlp Dener PetscFunctionBegin; 2869566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&auglag->subsolver)); 287661095bbSAlp Dener if (tao->setupcalled) { 2889566063dSJacob Faibussowitsch PetscCall(VecDestroy(&auglag->Xwork)); /* opt work */ 289661095bbSAlp Dener if (tao->eq_constrained) { 2909566063dSJacob Faibussowitsch PetscCall(VecDestroy(&auglag->Ye)); /* equality multipliers */ 2919566063dSJacob Faibussowitsch PetscCall(VecDestroy(&auglag->Cework)); /* equality work vector */ 292661095bbSAlp Dener } 293661095bbSAlp Dener if (tao->ineq_constrained) { 2949566063dSJacob Faibussowitsch PetscCall(VecDestroy(&auglag->Ps)); /* slack vars */ 295661095bbSAlp Dener auglag->Parr[0] = NULL; /* clear pointer to tao->solution, will be destroyed by TaoDestroy() shell */ 2969566063dSJacob Faibussowitsch PetscCall(PetscFree(auglag->Parr)); /* array of primal vectors */ 2979566063dSJacob Faibussowitsch PetscCall(VecDestroy(&auglag->LgradS)); /* slack grad */ 2989566063dSJacob Faibussowitsch PetscCall(VecDestroy(&auglag->Cizero)); /* zero vector for pointwise max */ 2999566063dSJacob Faibussowitsch PetscCall(VecDestroy(&auglag->Yi)); /* inequality multipliers */ 3009566063dSJacob Faibussowitsch PetscCall(VecDestroy(&auglag->Ciwork)); /* inequality work vector */ 3019566063dSJacob Faibussowitsch PetscCall(VecDestroy(&auglag->P)); /* combo primal */ 3029566063dSJacob Faibussowitsch PetscCall(ISDestroy(&auglag->Pis[0])); /* index set for X inside P */ 3039566063dSJacob Faibussowitsch PetscCall(ISDestroy(&auglag->Pis[1])); /* index set for S inside P */ 3049566063dSJacob Faibussowitsch PetscCall(PetscFree(auglag->Pis)); /* array of P index sets */ 3059566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&auglag->Pscatter[0])); 3069566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&auglag->Pscatter[1])); 3079566063dSJacob Faibussowitsch PetscCall(PetscFree(auglag->Pscatter)); 308661095bbSAlp Dener if (tao->eq_constrained) { 3099566063dSJacob Faibussowitsch PetscCall(VecDestroy(&auglag->Y)); /* combo multipliers */ 3109566063dSJacob Faibussowitsch PetscCall(PetscFree(auglag->Yarr)); /* array of dual vectors */ 3119566063dSJacob Faibussowitsch PetscCall(VecDestroy(&auglag->C)); /* combo constraints */ 3129566063dSJacob Faibussowitsch PetscCall(ISDestroy(&auglag->Yis[0])); /* index set for Ye inside Y */ 3139566063dSJacob Faibussowitsch PetscCall(ISDestroy(&auglag->Yis[1])); /* index set for Yi inside Y */ 3149566063dSJacob Faibussowitsch PetscCall(PetscFree(auglag->Yis)); 3159566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&auglag->Yscatter[0])); 3169566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&auglag->Yscatter[1])); 3179566063dSJacob Faibussowitsch PetscCall(PetscFree(auglag->Yscatter)); 318661095bbSAlp Dener } 319661095bbSAlp Dener } 320661095bbSAlp Dener if (tao->bounded) { 3219566063dSJacob Faibussowitsch PetscCall(VecDestroy(&auglag->PL)); /* lower bounds for subsolver */ 3229566063dSJacob Faibussowitsch PetscCall(VecDestroy(&auglag->PU)); /* upper bounds for subsolver */ 323661095bbSAlp Dener } 324661095bbSAlp Dener } 3252e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetType_C", NULL)); 3262e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMSetType_C", NULL)); 3272e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetSubsolver_C", NULL)); 3282e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMSetSubsolver_C", NULL)); 3292e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetMultipliers_C", NULL)); 3302e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMSetMultipliers_C", NULL)); 3312e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetPrimalIS_C", NULL)); 3322e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetDualIS_C", NULL)); 3339566063dSJacob Faibussowitsch PetscCall(PetscFree(tao->data)); 334661095bbSAlp Dener PetscFunctionReturn(0); 335661095bbSAlp Dener } 336661095bbSAlp Dener 337d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSetFromOptions_ALMM(Tao tao, PetscOptionItems *PetscOptionsObject) 338d71ae5a4SJacob Faibussowitsch { 339661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM *)tao->data; 340661095bbSAlp Dener PetscInt i; 341661095bbSAlp Dener 342661095bbSAlp Dener PetscFunctionBegin; 34371075aafSSatish Balay PetscOptionsHeadBegin(PetscOptionsObject, "Augmented Lagrangian multiplier method solves problems with general constraints by converting them into a sequence of unconstrained problems."); 3449566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_almm_mu_init", "initial penalty parameter", "", auglag->mu0, &auglag->mu0, NULL)); 3459566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_almm_mu_factor", "increase factor for the penalty parameter", "", auglag->mu_fac, &auglag->mu_fac, NULL)); 3469566063dSJacob Faibussowitsch 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)); 3479566063dSJacob Faibussowitsch 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)); 3489566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_almm_mu_max", "maximum safeguard for penalty parameter updates", "", auglag->mu_max, &auglag->mu_max, NULL)); 3499566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_almm_ye_min", "minimum safeguard for equality multiplier updates", "", auglag->ye_min, &auglag->ye_min, NULL)); 3509566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_almm_ye_max", "maximum safeguard for equality multipliers updates", "", auglag->ye_max, &auglag->ye_max, NULL)); 3519566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_almm_yi_min", "minimum safeguard for inequality multipliers updates", "", auglag->yi_min, &auglag->yi_min, NULL)); 3529566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_almm_yi_max", "maximum safeguard for inequality multipliers updates", "", auglag->yi_max, &auglag->yi_max, NULL)); 3539566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-tao_almm_type", "augmented Lagrangian formulation type for the subproblem", "TaoALMMType", TaoALMMTypes, (PetscEnum)auglag->type, (PetscEnum *)&auglag->type, NULL)); 354d0609cedSBarry Smith PetscOptionsHeadEnd(); 3559566063dSJacob Faibussowitsch PetscCall(TaoSetOptionsPrefix(auglag->subsolver, ((PetscObject)tao)->prefix)); 3569566063dSJacob Faibussowitsch PetscCall(TaoAppendOptionsPrefix(auglag->subsolver, "tao_almm_subsolver_")); 3579566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(auglag->subsolver)); 358661095bbSAlp Dener for (i = 0; i < tao->numbermonitors; i++) { 3599566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)tao->monitorcontext[i])); 3609566063dSJacob Faibussowitsch PetscCall(TaoSetMonitor(auglag->subsolver, tao->monitor[i], tao->monitorcontext[i], tao->monitordestroy[i])); 361ad540459SPierre Jolivet if (tao->monitor[i] == TaoMonitorDefault || tao->monitor[i] == TaoDefaultCMonitor || tao->monitor[i] == TaoDefaultGMonitor || tao->monitor[i] == TaoDefaultSMonitor) auglag->info = PETSC_TRUE; 362661095bbSAlp Dener } 363661095bbSAlp Dener PetscFunctionReturn(0); 364661095bbSAlp Dener } 365661095bbSAlp Dener 366661095bbSAlp Dener /* -------------------------------------------------------- */ 367661095bbSAlp Dener 368661095bbSAlp Dener /*MC 369*c0298470SBarry Smith TAOALMM - Augmented Lagrangian multiplier method for solving nonlinear optimization problems with general constraints. 370661095bbSAlp Dener 371661095bbSAlp Dener Options Database Keys: 372661095bbSAlp Dener + -tao_almm_mu_init <real> - initial penalty parameter (default: 10.) 373661095bbSAlp Dener . -tao_almm_mu_factor <real> - increase factor for the penalty parameter (default: 100.) 374661095bbSAlp Dener . -tao_almm_mu_max <real> - maximum safeguard for penalty parameter updates (default: 1.e20) 375661095bbSAlp Dener . -tao_almm_mu_power_good <real> - exponential for penalty parameter when multiplier update is accepted (default: 0.9) 376661095bbSAlp Dener . -tao_almm_mu_power_bad <real> - exponential for penalty parameter when multiplier update is rejected (default: 0.1) 377661095bbSAlp Dener . -tao_almm_ye_min <real> - minimum safeguard for equality multiplier updates (default: -1.e20) 378661095bbSAlp Dener . -tao_almm_ye_max <real> - maximum safeguard for equality multiplier updates (default: 1.e20) 379661095bbSAlp Dener . -tao_almm_yi_min <real> - minimum safeguard for inequality multiplier updates (default: -1.e20) 380661095bbSAlp Dener . -tao_almm_yi_max <real> - maximum safeguard for inequality multiplier updates (default: 1.e20) 381661095bbSAlp Dener - -tao_almm_type <classic,phr> - change formulation of the augmented Lagrangian merit function for the subproblem (default: classic) 382661095bbSAlp Dener 383661095bbSAlp Dener Level: beginner 384661095bbSAlp Dener 385661095bbSAlp Dener Notes: 386661095bbSAlp Dener This method converts a constrained problem into a sequence of unconstrained problems via the augmented 387661095bbSAlp Dener Lagrangian merit function. Bound constraints are pushed down to the subproblem without any modifications. 388661095bbSAlp Dener 389661095bbSAlp Dener Two formulations are offered for the subproblem: canonical Hestenes-Powell augmented Lagrangian with slack 390661095bbSAlp Dener variables for inequality constraints, and a slack-less Powell-Hestenes-Rockafellar (PHR) formulation utilizing a 391661095bbSAlp Dener pointwise max() penalty on inequality constraints. The canonical augmented Lagrangian formulation typically 392661095bbSAlp Dener converges faster for most problems. However, PHR may be desirable for problems featuring a large number 393661095bbSAlp Dener of inequality constraints because it avoids inflating the size of the subproblem with slack variables. 394661095bbSAlp Dener 395661095bbSAlp Dener The subproblem is solved using a nested first-order TAO solver. The user can retrieve a pointer to 396*c0298470SBarry Smith the subsolver via `TaoALMMGetSubsolver()` or pass command line arguments to it using the 397*c0298470SBarry Smith "-tao_almm_subsolver_" prefix. Currently, `TAOALMM` does not support second-order methods for the 398661095bbSAlp Dener subproblem. It is also highly recommended that the subsolver chosen by the user utilize a trust-region 399*c0298470SBarry Smith strategy for globalization (default: `TAOBQNKTR`) especially if the outer problem features bound constraints. 400661095bbSAlp Dener 401661095bbSAlp Dener .vb 402661095bbSAlp Dener while unconverged 403661095bbSAlp Dener solve argmin_x L(x) s.t. l <= x <= u 404661095bbSAlp Dener if ||c|| <= y_tol 405661095bbSAlp Dener if ||c|| <= c_tol && ||Lgrad|| <= g_tol: 406661095bbSAlp Dener problem converged, return solution 407661095bbSAlp Dener else 408661095bbSAlp Dener constraints sufficiently improved 409661095bbSAlp Dener update multipliers and tighten tolerances 410661095bbSAlp Dener endif 411661095bbSAlp Dener else 412661095bbSAlp Dener constraints did not improve 413661095bbSAlp Dener update penalty and loosen tolerances 414661095bbSAlp Dener endif 415661095bbSAlp Dener endwhile 416661095bbSAlp Dener .ve 417661095bbSAlp Dener 418*c0298470SBarry Smith .seealso: `TAOALMM`, `Tao`, `TaoALMMGetType()`, `TaoALMMSetType()`, `TaoALMMSetSubsolver()`, `TaoALMMGetSubsolver()`, 419db781477SPatrick Sanan `TaoALMMGetMultipliers()`, `TaoALMMSetMultipliers()`, `TaoALMMGetPrimalIS()`, `TaoALMMGetDualIS()` 420661095bbSAlp Dener M*/ 421d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TaoCreate_ALMM(Tao tao) 422d71ae5a4SJacob Faibussowitsch { 423661095bbSAlp Dener TAO_ALMM *auglag; 424661095bbSAlp Dener 425661095bbSAlp Dener PetscFunctionBegin; 4264dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&auglag)); 427661095bbSAlp Dener 428661095bbSAlp Dener tao->ops->destroy = TaoDestroy_ALMM; 429661095bbSAlp Dener tao->ops->setup = TaoSetUp_ALMM; 430661095bbSAlp Dener tao->ops->setfromoptions = TaoSetFromOptions_ALMM; 431661095bbSAlp Dener tao->ops->view = TaoView_ALMM; 432661095bbSAlp Dener tao->ops->solve = TaoSolve_ALMM; 433661095bbSAlp Dener 434661095bbSAlp Dener tao->gatol = 1.e-5; 435661095bbSAlp Dener tao->grtol = 0.0; 436661095bbSAlp Dener tao->gttol = 0.0; 437661095bbSAlp Dener tao->catol = 1.e-5; 438661095bbSAlp Dener tao->crtol = 0.0; 439661095bbSAlp Dener 440661095bbSAlp Dener tao->data = (void *)auglag; 441661095bbSAlp Dener auglag->parent = tao; 442661095bbSAlp Dener auglag->mu0 = 10.0; 443661095bbSAlp Dener auglag->mu = auglag->mu0; 444661095bbSAlp Dener auglag->mu_fac = 10.0; 445661095bbSAlp Dener auglag->mu_max = PETSC_INFINITY; 446661095bbSAlp Dener auglag->mu_pow_good = 0.9; 447661095bbSAlp Dener auglag->mu_pow_bad = 0.1; 448661095bbSAlp Dener auglag->ye_min = PETSC_NINFINITY; 449661095bbSAlp Dener auglag->ye_max = PETSC_INFINITY; 450661095bbSAlp Dener auglag->yi_min = PETSC_NINFINITY; 451661095bbSAlp Dener auglag->yi_max = PETSC_INFINITY; 452661095bbSAlp Dener auglag->ytol0 = 0.1 / PetscPowReal(auglag->mu0, auglag->mu_pow_bad); 453661095bbSAlp Dener auglag->ytol = auglag->ytol0; 454661095bbSAlp Dener auglag->gtol0 = 1.0 / auglag->mu0; 455661095bbSAlp Dener auglag->gtol = auglag->gtol0; 456661095bbSAlp Dener 457661095bbSAlp Dener auglag->sub_obj = TaoALMMComputeAugLagAndGradient_Private; 458661095bbSAlp Dener auglag->type = TAO_ALMM_CLASSIC; 459661095bbSAlp Dener auglag->info = PETSC_FALSE; 460661095bbSAlp Dener 4619566063dSJacob Faibussowitsch PetscCall(TaoCreate(PetscObjectComm((PetscObject)tao), &auglag->subsolver)); 4629566063dSJacob Faibussowitsch PetscCall(TaoSetType(auglag->subsolver, TAOBQNKTR)); 4639566063dSJacob Faibussowitsch PetscCall(TaoSetTolerances(auglag->subsolver, auglag->gtol, 0.0, 0.0)); 4649566063dSJacob Faibussowitsch PetscCall(TaoSetMaximumIterations(auglag->subsolver, 1000)); 4659566063dSJacob Faibussowitsch PetscCall(TaoSetMaximumFunctionEvaluations(auglag->subsolver, 10000)); 4669566063dSJacob Faibussowitsch PetscCall(TaoSetFunctionLowerBound(auglag->subsolver, PETSC_NINFINITY)); 4679566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)auglag->subsolver, (PetscObject)tao, 1)); 468661095bbSAlp Dener 4699566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetType_C", TaoALMMGetType_Private)); 4709566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMSetType_C", TaoALMMSetType_Private)); 4719566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetSubsolver_C", TaoALMMGetSubsolver_Private)); 4729566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMSetSubsolver_C", TaoALMMSetSubsolver_Private)); 4739566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetMultipliers_C", TaoALMMGetMultipliers_Private)); 4749566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMSetMultipliers_C", TaoALMMSetMultipliers_Private)); 4759566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetPrimalIS_C", TaoALMMGetPrimalIS_Private)); 4769566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetDualIS_C", TaoALMMGetDualIS_Private)); 477661095bbSAlp Dener PetscFunctionReturn(0); 478661095bbSAlp Dener } 479661095bbSAlp Dener 480d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoALMMCombinePrimal_Private(Tao tao, Vec X, Vec S, Vec P) 481d71ae5a4SJacob Faibussowitsch { 482661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM *)tao->data; 483661095bbSAlp Dener 484661095bbSAlp Dener PetscFunctionBegin; 485661095bbSAlp Dener if (tao->ineq_constrained) { 4869566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(auglag->Pscatter[0], X, P, INSERT_VALUES, SCATTER_REVERSE)); 4879566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(auglag->Pscatter[0], X, P, INSERT_VALUES, SCATTER_REVERSE)); 4889566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(auglag->Pscatter[1], S, P, INSERT_VALUES, SCATTER_REVERSE)); 4899566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(auglag->Pscatter[1], S, P, INSERT_VALUES, SCATTER_REVERSE)); 490661095bbSAlp Dener } else { 4919566063dSJacob Faibussowitsch PetscCall(VecCopy(X, P)); 492661095bbSAlp Dener } 493661095bbSAlp Dener PetscFunctionReturn(0); 494661095bbSAlp Dener } 495661095bbSAlp Dener 496d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoALMMCombineDual_Private(Tao tao, Vec EQ, Vec IN, Vec Y) 497d71ae5a4SJacob Faibussowitsch { 498661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM *)tao->data; 499661095bbSAlp Dener 500661095bbSAlp Dener PetscFunctionBegin; 501661095bbSAlp Dener if (tao->eq_constrained) { 502661095bbSAlp Dener if (tao->ineq_constrained) { 5039566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(auglag->Yscatter[0], EQ, Y, INSERT_VALUES, SCATTER_REVERSE)); 5049566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(auglag->Yscatter[0], EQ, Y, INSERT_VALUES, SCATTER_REVERSE)); 5059566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(auglag->Yscatter[1], IN, Y, INSERT_VALUES, SCATTER_REVERSE)); 5069566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(auglag->Yscatter[1], IN, Y, INSERT_VALUES, SCATTER_REVERSE)); 507661095bbSAlp Dener } else { 5089566063dSJacob Faibussowitsch PetscCall(VecCopy(EQ, Y)); 509661095bbSAlp Dener } 510661095bbSAlp Dener } else { 5119566063dSJacob Faibussowitsch PetscCall(VecCopy(IN, Y)); 512661095bbSAlp Dener } 513661095bbSAlp Dener PetscFunctionReturn(0); 514661095bbSAlp Dener } 515661095bbSAlp Dener 516d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoALMMSplitPrimal_Private(Tao tao, Vec P, Vec X, Vec S) 517d71ae5a4SJacob Faibussowitsch { 518661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM *)tao->data; 519661095bbSAlp Dener 520661095bbSAlp Dener PetscFunctionBegin; 521661095bbSAlp Dener if (tao->ineq_constrained) { 5229566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(auglag->Pscatter[0], P, X, INSERT_VALUES, SCATTER_FORWARD)); 5239566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(auglag->Pscatter[0], P, X, INSERT_VALUES, SCATTER_FORWARD)); 5249566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(auglag->Pscatter[1], P, S, INSERT_VALUES, SCATTER_FORWARD)); 5259566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(auglag->Pscatter[1], P, S, INSERT_VALUES, SCATTER_FORWARD)); 526661095bbSAlp Dener } else { 5279566063dSJacob Faibussowitsch PetscCall(VecCopy(P, X)); 528661095bbSAlp Dener } 529661095bbSAlp Dener PetscFunctionReturn(0); 530661095bbSAlp Dener } 531661095bbSAlp Dener 532661095bbSAlp Dener /* this assumes that the latest constraints are stored in Ce and Ci, and also combined in C */ 533d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoALMMComputeOptimalityNorms_Private(Tao tao) 534d71ae5a4SJacob Faibussowitsch { 535661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM *)tao->data; 536661095bbSAlp Dener 537661095bbSAlp Dener PetscFunctionBegin; 538661095bbSAlp Dener /* if bounded, project the gradient */ 5391baa6e33SBarry Smith if (tao->bounded) PetscCall(VecBoundGradientProjection(auglag->LgradX, auglag->Px, tao->XL, tao->XU, auglag->LgradX)); 540661095bbSAlp Dener if (auglag->type == TAO_ALMM_PHR) { 5419566063dSJacob Faibussowitsch PetscCall(VecNorm(auglag->LgradX, NORM_INFINITY, &auglag->gnorm)); 542661095bbSAlp Dener auglag->cenorm = 0.0; 54348a46eb9SPierre Jolivet if (tao->eq_constrained) PetscCall(VecNorm(auglag->Ce, NORM_INFINITY, &auglag->cenorm)); 544661095bbSAlp Dener auglag->cinorm = 0.0; 545661095bbSAlp Dener if (tao->ineq_constrained) { 5469566063dSJacob Faibussowitsch PetscCall(VecCopy(auglag->Yi, auglag->Ciwork)); 5479566063dSJacob Faibussowitsch PetscCall(VecScale(auglag->Ciwork, -1.0 / auglag->mu)); 5489566063dSJacob Faibussowitsch PetscCall(VecPointwiseMax(auglag->Ciwork, auglag->Ci, auglag->Ciwork)); 5499566063dSJacob Faibussowitsch PetscCall(VecNorm(auglag->Ciwork, NORM_INFINITY, &auglag->cinorm)); 550661095bbSAlp Dener } 551661095bbSAlp Dener auglag->cnorm_old = auglag->cnorm; 552661095bbSAlp Dener auglag->cnorm = PetscMax(auglag->cenorm, auglag->cinorm); 553661095bbSAlp Dener auglag->ytol = auglag->ytol0 * auglag->cnorm_old; 554661095bbSAlp Dener } else { 5559566063dSJacob Faibussowitsch PetscCall(VecNorm(auglag->LgradX, NORM_2, &auglag->gnorm)); 5569566063dSJacob Faibussowitsch PetscCall(VecNorm(auglag->C, NORM_2, &auglag->cnorm)); 557661095bbSAlp Dener } 558661095bbSAlp Dener PetscFunctionReturn(0); 559661095bbSAlp Dener } 560661095bbSAlp Dener 561d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoALMMEvaluateIterate_Private(Tao tao, Vec P) 562d71ae5a4SJacob Faibussowitsch { 563661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM *)tao->data; 564661095bbSAlp Dener 565661095bbSAlp Dener PetscFunctionBegin; 566661095bbSAlp Dener /* split solution into primal and slack components */ 5679566063dSJacob Faibussowitsch PetscCall(TaoALMMSplitPrimal_Private(tao, auglag->P, auglag->Px, auglag->Ps)); 568661095bbSAlp Dener 569661095bbSAlp Dener /* compute f, df/dx and the constraints */ 5709566063dSJacob Faibussowitsch PetscCall(TaoComputeObjectiveAndGradient(tao, auglag->Px, &auglag->fval, auglag->LgradX)); 571661095bbSAlp Dener if (tao->eq_constrained) { 5729566063dSJacob Faibussowitsch PetscCall(TaoComputeEqualityConstraints(tao, auglag->Px, auglag->Ce)); 5739566063dSJacob Faibussowitsch PetscCall(TaoComputeJacobianEquality(tao, auglag->Px, auglag->Ae, auglag->Ae)); 574661095bbSAlp Dener } 575661095bbSAlp Dener if (tao->ineq_constrained) { 5769566063dSJacob Faibussowitsch PetscCall(TaoComputeInequalityConstraints(tao, auglag->Px, auglag->Ci)); 5779566063dSJacob Faibussowitsch PetscCall(TaoComputeJacobianInequality(tao, auglag->Px, auglag->Ai, auglag->Ai)); 578661095bbSAlp Dener switch (auglag->type) { 5797721a36fSStefano Zampini case TAO_ALMM_CLASSIC: 580661095bbSAlp Dener /* classic formulation converts inequality to equality constraints via slack variables */ 5819566063dSJacob Faibussowitsch PetscCall(VecAXPY(auglag->Ci, -1.0, auglag->Ps)); 582661095bbSAlp Dener break; 5837721a36fSStefano Zampini case TAO_ALMM_PHR: 584661095bbSAlp Dener /* PHR is based on Ci <= 0 while TAO defines Ci >= 0 so we hit it with a negative sign */ 5859566063dSJacob Faibussowitsch PetscCall(VecScale(auglag->Ci, -1.0)); 5869566063dSJacob Faibussowitsch PetscCall(MatScale(auglag->Ai, -1.0)); 587661095bbSAlp Dener break; 588d71ae5a4SJacob Faibussowitsch default: 589d71ae5a4SJacob Faibussowitsch break; 590661095bbSAlp Dener } 591661095bbSAlp Dener } 592661095bbSAlp Dener /* combine constraints into one vector */ 5939566063dSJacob Faibussowitsch PetscCall(TaoALMMCombineDual_Private(tao, auglag->Ce, auglag->Ci, auglag->C)); 594661095bbSAlp Dener PetscFunctionReturn(0); 595661095bbSAlp Dener } 596661095bbSAlp Dener 597661095bbSAlp Dener /* 598661095bbSAlp Dener Lphr = f + 0.5*mu*[ (Ce + Ye/mu)^T (Ce + Ye/mu) + pmin(0, Ci + Yi/mu)^T pmin(0, Ci + Yi/mu)] 599661095bbSAlp Dener 600661095bbSAlp Dener dLphr/dX = dF/dX + mu*[ (Ce + Ye/mu)^T Ae + pmin(0, Ci + Yi/mu)^T Ai] 601661095bbSAlp Dener 602661095bbSAlp Dener dLphr/dS = 0 603661095bbSAlp Dener */ 604d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoALMMComputePHRLagAndGradient_Private(Tao tao) 605d71ae5a4SJacob Faibussowitsch { 606661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM *)tao->data; 607661095bbSAlp Dener PetscReal eq_norm = 0.0, ineq_norm = 0.0; 608661095bbSAlp Dener 609661095bbSAlp Dener PetscFunctionBegin; 6109566063dSJacob Faibussowitsch PetscCall(TaoALMMEvaluateIterate_Private(tao, auglag->P)); 611661095bbSAlp Dener if (tao->eq_constrained) { 612661095bbSAlp Dener /* Ce_work = mu*(Ce + Ye/mu) */ 6139566063dSJacob Faibussowitsch PetscCall(VecWAXPY(auglag->Cework, 1.0 / auglag->mu, auglag->Ye, auglag->Ce)); 6149566063dSJacob Faibussowitsch PetscCall(VecDot(auglag->Cework, auglag->Cework, &eq_norm)); /* contribution to scalar Lagrangian */ 6159566063dSJacob Faibussowitsch PetscCall(VecScale(auglag->Cework, auglag->mu)); 616661095bbSAlp Dener /* dL/dX += mu*(Ce + Ye/mu)^T Ae */ 6179566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(auglag->Ae, auglag->Cework, auglag->LgradX, auglag->LgradX)); 618661095bbSAlp Dener } 619661095bbSAlp Dener if (tao->ineq_constrained) { 620661095bbSAlp Dener /* Ci_work = mu * pmax(0, Ci + Yi/mu) where pmax() is pointwise max() */ 6219566063dSJacob Faibussowitsch PetscCall(VecWAXPY(auglag->Ciwork, 1.0 / auglag->mu, auglag->Yi, auglag->Ci)); 6229566063dSJacob Faibussowitsch PetscCall(VecPointwiseMax(auglag->Ciwork, auglag->Cizero, auglag->Ciwork)); 6239566063dSJacob Faibussowitsch PetscCall(VecDot(auglag->Ciwork, auglag->Ciwork, &ineq_norm)); /* contribution to scalar Lagrangian */ 624661095bbSAlp Dener /* dL/dX += mu * pmax(0, Ci + Yi/mu)^T Ai */ 6259566063dSJacob Faibussowitsch PetscCall(VecScale(auglag->Ciwork, auglag->mu)); 6269566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(auglag->Ai, auglag->Ciwork, auglag->LgradX, auglag->LgradX)); 627a5b23f4aSJose E. Roman /* dL/dS = 0 because there are no slacks in PHR */ 6289566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(auglag->LgradS)); 629661095bbSAlp Dener } 630661095bbSAlp Dener /* combine gradient together */ 6319566063dSJacob Faibussowitsch PetscCall(TaoALMMCombinePrimal_Private(tao, auglag->LgradX, auglag->LgradS, auglag->G)); 632661095bbSAlp Dener /* compute L = f + 0.5 * mu * [(Ce + Ye/mu)^T (Ce + Ye/mu) + pmax(0, Ci + Yi/mu)^T pmax(0, Ci + Yi/mu)] */ 6335f80ce2aSJacob Faibussowitsch auglag->Lval = auglag->fval + 0.5 * auglag->mu * (eq_norm + ineq_norm); 634661095bbSAlp Dener PetscFunctionReturn(0); 635661095bbSAlp Dener } 636661095bbSAlp Dener 637661095bbSAlp Dener /* 638661095bbSAlp Dener Lc = F + Ye^TCe + Yi^T(Ci - S) + 0.5*mu*[Ce^TCe + (Ci - S)^T(Ci - S)] 639661095bbSAlp Dener 640661095bbSAlp Dener dLc/dX = dF/dX + Ye^TAe + Yi^TAi + 0.5*mu*[Ce^TAe + (Ci - S)^TAi] 641661095bbSAlp Dener 642661095bbSAlp Dener dLc/dS = -[Yi + mu*(Ci - S)] 643661095bbSAlp Dener */ 644d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoALMMComputeAugLagAndGradient_Private(Tao tao) 645d71ae5a4SJacob Faibussowitsch { 646661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM *)tao->data; 647661095bbSAlp Dener PetscReal yeTce = 0.0, yiTcims = 0.0, ceTce = 0.0, cimsTcims = 0.0; 648661095bbSAlp Dener 649661095bbSAlp Dener PetscFunctionBegin; 6509566063dSJacob Faibussowitsch PetscCall(TaoALMMEvaluateIterate_Private(tao, auglag->P)); 651661095bbSAlp Dener if (tao->eq_constrained) { 652661095bbSAlp Dener /* compute scalar contributions */ 6539566063dSJacob Faibussowitsch PetscCall(VecDot(auglag->Ye, auglag->Ce, &yeTce)); 6549566063dSJacob Faibussowitsch PetscCall(VecDot(auglag->Ce, auglag->Ce, &ceTce)); 655661095bbSAlp Dener /* dL/dX += ye^T Ae */ 6569566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(auglag->Ae, auglag->Ye, auglag->LgradX, auglag->LgradX)); 657661095bbSAlp Dener /* dL/dX += mu * ce^T Ae */ 6589566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(auglag->Ae, auglag->Ce, auglag->Xwork)); 6599566063dSJacob Faibussowitsch PetscCall(VecAXPY(auglag->LgradX, auglag->mu, auglag->Xwork)); 660661095bbSAlp Dener } 661661095bbSAlp Dener if (tao->ineq_constrained) { 662661095bbSAlp Dener /* compute scalar contributions */ 6639566063dSJacob Faibussowitsch PetscCall(VecDot(auglag->Yi, auglag->Ci, &yiTcims)); 6649566063dSJacob Faibussowitsch PetscCall(VecDot(auglag->Ci, auglag->Ci, &cimsTcims)); 665661095bbSAlp Dener /* dL/dX += yi^T Ai */ 6669566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(auglag->Ai, auglag->Yi, auglag->LgradX, auglag->LgradX)); 667661095bbSAlp Dener /* dL/dX += mu * (ci - s)^T Ai */ 6689566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(auglag->Ai, auglag->Ci, auglag->Xwork)); 6699566063dSJacob Faibussowitsch PetscCall(VecAXPY(auglag->LgradX, auglag->mu, auglag->Xwork)); 670661095bbSAlp Dener /* dL/dS = -[yi + mu*(ci - s)] */ 6719566063dSJacob Faibussowitsch PetscCall(VecWAXPY(auglag->LgradS, auglag->mu, auglag->Ci, auglag->Yi)); 6729566063dSJacob Faibussowitsch PetscCall(VecScale(auglag->LgradS, -1.0)); 673661095bbSAlp Dener } 674661095bbSAlp Dener /* combine gradient together */ 6759566063dSJacob Faibussowitsch PetscCall(TaoALMMCombinePrimal_Private(tao, auglag->LgradX, auglag->LgradS, auglag->G)); 676661095bbSAlp Dener /* compute L = f + ye^T ce + yi^T (ci - s) + 0.5*mu*||ce||^2 + 0.5*mu*||ci - s||^2 */ 677661095bbSAlp Dener auglag->Lval = auglag->fval + yeTce + yiTcims + 0.5 * auglag->mu * (ceTce + cimsTcims); 678661095bbSAlp Dener PetscFunctionReturn(0); 679661095bbSAlp Dener } 680661095bbSAlp Dener 681d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoALMMSubsolverObjective_Private(Tao tao, Vec P, PetscReal *Lval, void *ctx) 682d71ae5a4SJacob Faibussowitsch { 6837721a36fSStefano Zampini TAO_ALMM *auglag = (TAO_ALMM *)ctx; 6847721a36fSStefano Zampini 6857721a36fSStefano Zampini PetscFunctionBegin; 6869566063dSJacob Faibussowitsch PetscCall(VecCopy(P, auglag->P)); 6879566063dSJacob Faibussowitsch PetscCall((*auglag->sub_obj)(auglag->parent)); 6887721a36fSStefano Zampini *Lval = auglag->Lval; 6897721a36fSStefano Zampini PetscFunctionReturn(0); 6907721a36fSStefano Zampini } 6917721a36fSStefano Zampini 692d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoALMMSubsolverObjectiveAndGradient_Private(Tao tao, Vec P, PetscReal *Lval, Vec G, void *ctx) 693d71ae5a4SJacob Faibussowitsch { 694661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM *)ctx; 695661095bbSAlp Dener 696661095bbSAlp Dener PetscFunctionBegin; 6979566063dSJacob Faibussowitsch PetscCall(VecCopy(P, auglag->P)); 6989566063dSJacob Faibussowitsch PetscCall((*auglag->sub_obj)(auglag->parent)); 6999566063dSJacob Faibussowitsch PetscCall(VecCopy(auglag->G, G)); 700661095bbSAlp Dener *Lval = auglag->Lval; 701661095bbSAlp Dener PetscFunctionReturn(0); 702661095bbSAlp Dener } 703