1954e39ddSJose E. Roman #include <../src/tao/constrained/impls/almm/almm.h> /*I "petsctao.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)); 259a09327dSAlp Dener PetscCall(VecSet(auglag->Yi, 1.0)); 26661095bbSAlp Dener } 279a09327dSAlp Dener if (tao->eq_constrained) PetscCall(VecSet(auglag->Ye, 1.0)); 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 */ 54*08fc2fa2SPaul T. Kühner if (PetscAbsReal(auglag->cenorm + auglag->cinorm) == 0.0) auglag->mu = 10.0; 55*08fc2fa2SPaul T. Kühner else auglag->mu = PetscMax(1.e-6, PetscMin(10.0, 2.0 * PetscAbsReal(auglag->fval) / (auglag->cenorm + auglag->cinorm))); 56661095bbSAlp Dener break; 57d71ae5a4SJacob Faibussowitsch default: 58d71ae5a4SJacob Faibussowitsch break; 59661095bbSAlp Dener } 60661095bbSAlp Dener auglag->gtol = auglag->gtol0; 6163a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(tao, "Initial penalty: %.2g\n", (double)auglag->mu)); 62661095bbSAlp Dener 63661095bbSAlp Dener /* start aug-lag outer loop */ 64661095bbSAlp Dener while (tao->reason == TAO_CONTINUE_ITERATING) { 65661095bbSAlp Dener ++tao->niter; 66661095bbSAlp Dener /* update subsolver tolerance */ 6763a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(tao, "Subsolver tolerance: ||G|| <= %e\n", (double)auglag->gtol)); 689566063dSJacob Faibussowitsch PetscCall(TaoSetTolerances(auglag->subsolver, auglag->gtol, 0.0, 0.0)); 69661095bbSAlp Dener /* solve the bound-constrained or unconstrained subproblem */ 709a09327dSAlp Dener PetscCall(VecCopy(auglag->P, auglag->subsolver->solution)); 719566063dSJacob Faibussowitsch PetscCall(TaoSolve(auglag->subsolver)); 729a09327dSAlp Dener PetscCall(VecCopy(auglag->subsolver->solution, auglag->P)); 739566063dSJacob Faibussowitsch PetscCall(TaoGetConvergedReason(auglag->subsolver, &reason)); 74661095bbSAlp Dener tao->ksp_its += auglag->subsolver->ksp_its; 7548a46eb9SPierre Jolivet if (reason != TAO_CONVERGED_GATOL) PetscCall(PetscInfo(tao, "Subsolver failed to converge, reason: %s\n", TaoConvergedReasons[reason])); 76661095bbSAlp Dener /* evaluate solution and test convergence */ 779566063dSJacob Faibussowitsch PetscCall((*auglag->sub_obj)(tao)); 789566063dSJacob Faibussowitsch PetscCall(TaoALMMComputeOptimalityNorms_Private(tao)); 79661095bbSAlp Dener /* decide whether to update multipliers or not */ 80661095bbSAlp Dener updated = 0.0; 81661095bbSAlp Dener if (auglag->cnorm <= auglag->ytol) { 8263a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(tao, "Multipliers updated: ||C|| <= %e\n", (double)auglag->ytol)); 83661095bbSAlp Dener /* constraints are good, update multipliers and convergence tolerances */ 84661095bbSAlp Dener if (tao->eq_constrained) { 859566063dSJacob Faibussowitsch PetscCall(VecAXPY(auglag->Ye, auglag->mu, auglag->Ce)); 869566063dSJacob Faibussowitsch PetscCall(VecSet(auglag->Cework, auglag->ye_max)); 879566063dSJacob Faibussowitsch PetscCall(VecPointwiseMin(auglag->Ye, auglag->Cework, auglag->Ye)); 889566063dSJacob Faibussowitsch PetscCall(VecSet(auglag->Cework, auglag->ye_min)); 899566063dSJacob Faibussowitsch PetscCall(VecPointwiseMax(auglag->Ye, auglag->Cework, auglag->Ye)); 90661095bbSAlp Dener } 91661095bbSAlp Dener if (tao->ineq_constrained) { 929566063dSJacob Faibussowitsch PetscCall(VecAXPY(auglag->Yi, auglag->mu, auglag->Ci)); 939566063dSJacob Faibussowitsch PetscCall(VecSet(auglag->Ciwork, auglag->yi_max)); 949566063dSJacob Faibussowitsch PetscCall(VecPointwiseMin(auglag->Yi, auglag->Ciwork, auglag->Yi)); 959566063dSJacob Faibussowitsch PetscCall(VecSet(auglag->Ciwork, auglag->yi_min)); 969566063dSJacob Faibussowitsch PetscCall(VecPointwiseMax(auglag->Yi, auglag->Ciwork, auglag->Yi)); 97661095bbSAlp Dener } 98661095bbSAlp Dener /* tolerances are updated only for non-PHR methods */ 99661095bbSAlp Dener if (auglag->type != TAO_ALMM_PHR) { 100661095bbSAlp Dener auglag->ytol = PetscMax(tao->catol, auglag->ytol / PetscPowReal(auglag->mu, auglag->mu_pow_good)); 101661095bbSAlp Dener auglag->gtol = PetscMax(tao->gatol, auglag->gtol / auglag->mu); 102661095bbSAlp Dener } 103661095bbSAlp Dener updated = 1.0; 104661095bbSAlp Dener } else { 105661095bbSAlp Dener /* constraints are bad, update penalty factor */ 106661095bbSAlp Dener auglag->mu = PetscMin(auglag->mu_max, auglag->mu_fac * auglag->mu); 107661095bbSAlp Dener /* tolerances are reset only for non-PHR methods */ 108661095bbSAlp Dener if (auglag->type != TAO_ALMM_PHR) { 109661095bbSAlp Dener auglag->ytol = PetscMax(tao->catol, 0.1 / PetscPowReal(auglag->mu, auglag->mu_pow_bad)); 110661095bbSAlp Dener auglag->gtol = PetscMax(tao->gatol, 1.0 / auglag->mu); 111661095bbSAlp Dener } 11263a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(tao, "Penalty increased: mu = %.2g\n", (double)auglag->mu)); 113661095bbSAlp Dener } 1149566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao, auglag->fval, auglag->gnorm, auglag->cnorm, tao->ksp_its)); 1159566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao, tao->niter, auglag->fval, auglag->gnorm, auglag->cnorm, updated)); 116dbbe0bcdSBarry Smith PetscUseTypeMethod(tao, convergencetest, tao->cnvP); 117661095bbSAlp Dener } 1183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 119661095bbSAlp Dener } 120661095bbSAlp Dener 121d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoView_ALMM(Tao tao, PetscViewer viewer) 122d71ae5a4SJacob Faibussowitsch { 123661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM *)tao->data; 124661095bbSAlp Dener PetscBool isascii; 125661095bbSAlp Dener 126661095bbSAlp Dener PetscFunctionBegin; 1279566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 1289566063dSJacob Faibussowitsch PetscCall(TaoView(auglag->subsolver, viewer)); 1299566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 1309566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 131661095bbSAlp Dener if (isascii) { 1329566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 1339566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "ALMM Formulation Type: %s\n", TaoALMMTypes[auglag->type])); 1349566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 135661095bbSAlp Dener } 1363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 137661095bbSAlp Dener } 138661095bbSAlp Dener 139d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoSetUp_ALMM(Tao tao) 140d71ae5a4SJacob Faibussowitsch { 141661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM *)tao->data; 142661095bbSAlp Dener VecType vec_type; 143661095bbSAlp Dener Vec SL, SU; 144661095bbSAlp Dener PetscBool is_cg = PETSC_FALSE, is_lmvm = PETSC_FALSE; 145661095bbSAlp Dener 146661095bbSAlp Dener PetscFunctionBegin; 14769d47153SPierre Jolivet PetscCheck(!tao->ineq_doublesided, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "TAOALMM does not support double-sided inequality constraint definition. Please restructure your inequality constraint to fit the form c(x) >= 0."); 1483c859ba3SBarry 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."); 1499566063dSJacob Faibussowitsch PetscCall(TaoComputeVariableBounds(tao)); 150661095bbSAlp Dener /* alias base vectors and create extras */ 1519566063dSJacob Faibussowitsch PetscCall(VecGetType(tao->solution, &vec_type)); 152661095bbSAlp Dener auglag->Px = tao->solution; 153661095bbSAlp Dener if (!tao->gradient) { /* base gradient */ 1549566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &tao->gradient)); 155661095bbSAlp Dener } 156661095bbSAlp Dener auglag->LgradX = tao->gradient; 157661095bbSAlp Dener if (!auglag->Xwork) { /* opt var work vector */ 1589566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &auglag->Xwork)); 159661095bbSAlp Dener } 160661095bbSAlp Dener if (tao->eq_constrained) { 161661095bbSAlp Dener auglag->Ce = tao->constraints_equality; 162661095bbSAlp Dener auglag->Ae = tao->jacobian_equality; 163661095bbSAlp Dener if (!auglag->Ye) { /* equality multipliers */ 1649566063dSJacob Faibussowitsch PetscCall(VecDuplicate(auglag->Ce, &auglag->Ye)); 165661095bbSAlp Dener } 16648a46eb9SPierre Jolivet if (!auglag->Cework) PetscCall(VecDuplicate(auglag->Ce, &auglag->Cework)); 167661095bbSAlp Dener } 168661095bbSAlp Dener if (tao->ineq_constrained) { 169661095bbSAlp Dener auglag->Ci = tao->constraints_inequality; 170661095bbSAlp Dener auglag->Ai = tao->jacobian_inequality; 171661095bbSAlp Dener if (!auglag->Yi) { /* inequality multipliers */ 1729566063dSJacob Faibussowitsch PetscCall(VecDuplicate(auglag->Ci, &auglag->Yi)); 173661095bbSAlp Dener } 17448a46eb9SPierre Jolivet if (!auglag->Ciwork) PetscCall(VecDuplicate(auglag->Ci, &auglag->Ciwork)); 175661095bbSAlp Dener if (!auglag->Cizero) { 1769566063dSJacob Faibussowitsch PetscCall(VecDuplicate(auglag->Ci, &auglag->Cizero)); 1779566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(auglag->Cizero)); 178661095bbSAlp Dener } 179661095bbSAlp Dener if (!auglag->Ps) { /* slack vars */ 1809566063dSJacob Faibussowitsch PetscCall(VecDuplicate(auglag->Ci, &auglag->Ps)); 181661095bbSAlp Dener } 182661095bbSAlp Dener if (!auglag->LgradS) { /* slack component of Lagrangian gradient */ 1839566063dSJacob Faibussowitsch PetscCall(VecDuplicate(auglag->Ci, &auglag->LgradS)); 184661095bbSAlp Dener } 185661095bbSAlp Dener /* create vector for combined primal space and the associated communication objects */ 186661095bbSAlp Dener if (!auglag->P) { 1879566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2, &auglag->Parr)); 1889371c9d4SSatish Balay auglag->Parr[0] = auglag->Px; 1899371c9d4SSatish Balay auglag->Parr[1] = auglag->Ps; 1909566063dSJacob Faibussowitsch PetscCall(VecConcatenate(2, auglag->Parr, &auglag->P, &auglag->Pis)); 1919566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2, &auglag->Pscatter)); 1929566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(auglag->P, auglag->Pis[0], auglag->Px, NULL, &auglag->Pscatter[0])); 1939566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(auglag->P, auglag->Pis[1], auglag->Ps, NULL, &auglag->Pscatter[1])); 194661095bbSAlp Dener } 195661095bbSAlp Dener if (tao->eq_constrained) { 196661095bbSAlp Dener /* create vector for combined dual space and the associated communication objects */ 197661095bbSAlp Dener if (!auglag->Y) { 1989566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2, &auglag->Yarr)); 1999371c9d4SSatish Balay auglag->Yarr[0] = auglag->Ye; 2009371c9d4SSatish Balay auglag->Yarr[1] = auglag->Yi; 2019566063dSJacob Faibussowitsch PetscCall(VecConcatenate(2, auglag->Yarr, &auglag->Y, &auglag->Yis)); 2029566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2, &auglag->Yscatter)); 2039566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(auglag->Y, auglag->Yis[0], auglag->Ye, NULL, &auglag->Yscatter[0])); 2049566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(auglag->Y, auglag->Yis[1], auglag->Yi, NULL, &auglag->Yscatter[1])); 205661095bbSAlp Dener } 20648a46eb9SPierre Jolivet if (!auglag->C) PetscCall(VecDuplicate(auglag->Y, &auglag->C)); 2079371c9d4SSatish Balay } else { 208ad540459SPierre Jolivet if (!auglag->C) auglag->C = auglag->Ci; 209ad540459SPierre Jolivet if (!auglag->Y) auglag->Y = auglag->Yi; 210661095bbSAlp Dener } 211661095bbSAlp Dener } else { 212ad540459SPierre Jolivet if (!auglag->P) auglag->P = auglag->Px; 213ad540459SPierre Jolivet if (!auglag->G) auglag->G = auglag->LgradX; 214ad540459SPierre Jolivet if (!auglag->C) auglag->C = auglag->Ce; 215ad540459SPierre Jolivet if (!auglag->Y) auglag->Y = auglag->Ye; 216661095bbSAlp Dener } 2179a09327dSAlp Dener /* create tao primal solution and gradient to interface with subsolver */ 2189a09327dSAlp Dener PetscCall(VecDuplicate(auglag->P, &auglag->Psub)); 2199a09327dSAlp Dener PetscCall(VecDuplicate(auglag->P, &auglag->G)); 220661095bbSAlp Dener /* initialize parameters */ 221661095bbSAlp Dener if (auglag->type == TAO_ALMM_PHR) { 222661095bbSAlp Dener auglag->mu_fac = 10.0; 223661095bbSAlp Dener auglag->yi_min = 0.0; 224661095bbSAlp Dener auglag->ytol0 = 0.5; 225661095bbSAlp Dener auglag->gtol0 = tao->gatol; 226606f75f6SBarry Smith if (tao->gatol != tao->default_gatol && tao->catol != tao->default_catol) { 2279566063dSJacob Faibussowitsch PetscCall(PetscInfo(tao, "TAOALMM with PHR: different gradient and constraint tolerances are not supported, setting catol = gatol\n")); 228661095bbSAlp Dener tao->catol = tao->gatol; 229661095bbSAlp Dener } 230661095bbSAlp Dener } 231661095bbSAlp Dener /* set the Lagrangian formulation type for the subsolver */ 232661095bbSAlp Dener switch (auglag->type) { 233d71ae5a4SJacob Faibussowitsch case TAO_ALMM_CLASSIC: 234d71ae5a4SJacob Faibussowitsch auglag->sub_obj = TaoALMMComputeAugLagAndGradient_Private; 235d71ae5a4SJacob Faibussowitsch break; 236d71ae5a4SJacob Faibussowitsch case TAO_ALMM_PHR: 237d71ae5a4SJacob Faibussowitsch auglag->sub_obj = TaoALMMComputePHRLagAndGradient_Private; 238d71ae5a4SJacob Faibussowitsch break; 239d71ae5a4SJacob Faibussowitsch default: 240d71ae5a4SJacob Faibussowitsch break; 241661095bbSAlp Dener } 242661095bbSAlp Dener /* set up the subsolver */ 2439a09327dSAlp Dener PetscCall(TaoSetSolution(auglag->subsolver, auglag->Psub)); 2449566063dSJacob Faibussowitsch PetscCall(TaoSetObjective(auglag->subsolver, TaoALMMSubsolverObjective_Private, (void *)auglag)); 2459566063dSJacob Faibussowitsch PetscCall(TaoSetObjectiveAndGradient(auglag->subsolver, NULL, TaoALMMSubsolverObjectiveAndGradient_Private, (void *)auglag)); 246661095bbSAlp Dener if (tao->bounded) { 247661095bbSAlp Dener /* make sure that the subsolver is a bound-constrained method */ 2489566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)auglag->subsolver, TAOCG, &is_cg)); 2499566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)auglag->subsolver, TAOLMVM, &is_lmvm)); 250661095bbSAlp Dener if (is_cg) { 2519566063dSJacob Faibussowitsch PetscCall(TaoSetType(auglag->subsolver, TAOBNCG)); 2529d3446b2SPierre Jolivet PetscCall(PetscInfo(tao, "TAOCG detected for bound-constrained problem, switching to TAOBNCG.\n")); 253661095bbSAlp Dener } 254661095bbSAlp Dener if (is_lmvm) { 2559566063dSJacob Faibussowitsch PetscCall(TaoSetType(auglag->subsolver, TAOBQNLS)); 2569d3446b2SPierre Jolivet PetscCall(PetscInfo(tao, "TAOLMVM detected for bound-constrained problem, switching to TAOBQNLS.\n")); 257661095bbSAlp Dener } 258661095bbSAlp Dener /* create lower and upper bound clone vectors for subsolver */ 25948a46eb9SPierre Jolivet if (!auglag->PL) PetscCall(VecDuplicate(auglag->P, &auglag->PL)); 26048a46eb9SPierre Jolivet if (!auglag->PU) PetscCall(VecDuplicate(auglag->P, &auglag->PU)); 261661095bbSAlp Dener if (tao->ineq_constrained) { 262661095bbSAlp Dener /* create lower and upper bounds for slack, set lower to 0 */ 2639566063dSJacob Faibussowitsch PetscCall(VecDuplicate(auglag->Ci, &SL)); 2649566063dSJacob Faibussowitsch PetscCall(VecSet(SL, 0.0)); 2659566063dSJacob Faibussowitsch PetscCall(VecDuplicate(auglag->Ci, &SU)); 2669566063dSJacob Faibussowitsch PetscCall(VecSet(SU, PETSC_INFINITY)); 267661095bbSAlp Dener /* combine opt var bounds with slack bounds */ 2689566063dSJacob Faibussowitsch PetscCall(TaoALMMCombinePrimal_Private(tao, tao->XL, SL, auglag->PL)); 2699566063dSJacob Faibussowitsch PetscCall(TaoALMMCombinePrimal_Private(tao, tao->XU, SU, auglag->PU)); 270661095bbSAlp Dener /* destroy work vectors */ 2719566063dSJacob Faibussowitsch PetscCall(VecDestroy(&SL)); 2729566063dSJacob Faibussowitsch PetscCall(VecDestroy(&SU)); 273661095bbSAlp Dener } else { 274661095bbSAlp Dener /* no inequality constraints, just copy bounds into the subsolver */ 2759566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->XL, auglag->PL)); 2769566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->XU, auglag->PU)); 277661095bbSAlp Dener } 2789566063dSJacob Faibussowitsch PetscCall(TaoSetVariableBounds(auglag->subsolver, auglag->PL, auglag->PU)); 279661095bbSAlp Dener } 2809566063dSJacob Faibussowitsch PetscCall(TaoSetUp(auglag->subsolver)); 2813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 282661095bbSAlp Dener } 283661095bbSAlp Dener 284d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoDestroy_ALMM(Tao tao) 285d71ae5a4SJacob Faibussowitsch { 286661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM *)tao->data; 287661095bbSAlp Dener 288661095bbSAlp Dener PetscFunctionBegin; 2899566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&auglag->subsolver)); 2909a09327dSAlp Dener PetscCall(VecDestroy(&auglag->Psub)); 2919a09327dSAlp Dener PetscCall(VecDestroy(&auglag->G)); 292661095bbSAlp Dener if (tao->setupcalled) { 2939a09327dSAlp Dener PetscCall(VecDestroy(&auglag->Xwork)); 294661095bbSAlp Dener if (tao->eq_constrained) { 2959566063dSJacob Faibussowitsch PetscCall(VecDestroy(&auglag->Ye)); /* equality multipliers */ 2969566063dSJacob Faibussowitsch PetscCall(VecDestroy(&auglag->Cework)); /* equality work vector */ 297661095bbSAlp Dener } 298661095bbSAlp Dener if (tao->ineq_constrained) { 2999566063dSJacob Faibussowitsch PetscCall(VecDestroy(&auglag->Ps)); /* slack vars */ 3009566063dSJacob Faibussowitsch PetscCall(PetscFree(auglag->Parr)); /* array of primal vectors */ 3019566063dSJacob Faibussowitsch PetscCall(VecDestroy(&auglag->LgradS)); /* slack grad */ 3029566063dSJacob Faibussowitsch PetscCall(VecDestroy(&auglag->Cizero)); /* zero vector for pointwise max */ 3039566063dSJacob Faibussowitsch PetscCall(VecDestroy(&auglag->Yi)); /* inequality multipliers */ 3049566063dSJacob Faibussowitsch PetscCall(VecDestroy(&auglag->Ciwork)); /* inequality work vector */ 3059a09327dSAlp Dener PetscCall(VecDestroy(&auglag->P)); /* combo primal solution */ 3069566063dSJacob Faibussowitsch PetscCall(ISDestroy(&auglag->Pis[0])); /* index set for X inside P */ 3079566063dSJacob Faibussowitsch PetscCall(ISDestroy(&auglag->Pis[1])); /* index set for S inside P */ 3089566063dSJacob Faibussowitsch PetscCall(PetscFree(auglag->Pis)); /* array of P index sets */ 3099566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&auglag->Pscatter[0])); 3109566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&auglag->Pscatter[1])); 3119566063dSJacob Faibussowitsch PetscCall(PetscFree(auglag->Pscatter)); 312661095bbSAlp Dener if (tao->eq_constrained) { 3139566063dSJacob Faibussowitsch PetscCall(VecDestroy(&auglag->Y)); /* combo multipliers */ 3149566063dSJacob Faibussowitsch PetscCall(PetscFree(auglag->Yarr)); /* array of dual vectors */ 3159566063dSJacob Faibussowitsch PetscCall(VecDestroy(&auglag->C)); /* combo constraints */ 3169566063dSJacob Faibussowitsch PetscCall(ISDestroy(&auglag->Yis[0])); /* index set for Ye inside Y */ 3179566063dSJacob Faibussowitsch PetscCall(ISDestroy(&auglag->Yis[1])); /* index set for Yi inside Y */ 3189566063dSJacob Faibussowitsch PetscCall(PetscFree(auglag->Yis)); 3199566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&auglag->Yscatter[0])); 3209566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&auglag->Yscatter[1])); 3219566063dSJacob Faibussowitsch PetscCall(PetscFree(auglag->Yscatter)); 322661095bbSAlp Dener } 323661095bbSAlp Dener } 324661095bbSAlp Dener if (tao->bounded) { 3259566063dSJacob Faibussowitsch PetscCall(VecDestroy(&auglag->PL)); /* lower bounds for subsolver */ 3269566063dSJacob Faibussowitsch PetscCall(VecDestroy(&auglag->PU)); /* upper bounds for subsolver */ 327661095bbSAlp Dener } 328661095bbSAlp Dener } 3292e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetType_C", NULL)); 3302e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMSetType_C", NULL)); 3312e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetSubsolver_C", NULL)); 3322e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMSetSubsolver_C", NULL)); 3332e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetMultipliers_C", NULL)); 3342e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMSetMultipliers_C", NULL)); 3352e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetPrimalIS_C", NULL)); 3362e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetDualIS_C", NULL)); 3379566063dSJacob Faibussowitsch PetscCall(PetscFree(tao->data)); 3383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 339661095bbSAlp Dener } 340661095bbSAlp Dener 341ce78bad3SBarry Smith static PetscErrorCode TaoSetFromOptions_ALMM(Tao tao, PetscOptionItems PetscOptionsObject) 342d71ae5a4SJacob Faibussowitsch { 343661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM *)tao->data; 344661095bbSAlp Dener PetscInt i; 345661095bbSAlp Dener 346661095bbSAlp Dener PetscFunctionBegin; 34771075aafSSatish Balay PetscOptionsHeadBegin(PetscOptionsObject, "Augmented Lagrangian multiplier method solves problems with general constraints by converting them into a sequence of unconstrained problems."); 3489566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_almm_mu_init", "initial penalty parameter", "", auglag->mu0, &auglag->mu0, NULL)); 3499566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_almm_mu_factor", "increase factor for the penalty parameter", "", auglag->mu_fac, &auglag->mu_fac, NULL)); 3509566063dSJacob 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)); 3519566063dSJacob 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)); 3529566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_almm_mu_max", "maximum safeguard for penalty parameter updates", "", auglag->mu_max, &auglag->mu_max, NULL)); 3539566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_almm_ye_min", "minimum safeguard for equality multiplier updates", "", auglag->ye_min, &auglag->ye_min, NULL)); 3549566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_almm_ye_max", "maximum safeguard for equality multipliers updates", "", auglag->ye_max, &auglag->ye_max, NULL)); 3559566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_almm_yi_min", "minimum safeguard for inequality multipliers updates", "", auglag->yi_min, &auglag->yi_min, NULL)); 3569566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_almm_yi_max", "maximum safeguard for inequality multipliers updates", "", auglag->yi_max, &auglag->yi_max, NULL)); 3579566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-tao_almm_type", "augmented Lagrangian formulation type for the subproblem", "TaoALMMType", TaoALMMTypes, (PetscEnum)auglag->type, (PetscEnum *)&auglag->type, NULL)); 358d0609cedSBarry Smith PetscOptionsHeadEnd(); 3599566063dSJacob Faibussowitsch PetscCall(TaoSetOptionsPrefix(auglag->subsolver, ((PetscObject)tao)->prefix)); 3609566063dSJacob Faibussowitsch PetscCall(TaoAppendOptionsPrefix(auglag->subsolver, "tao_almm_subsolver_")); 3619566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(auglag->subsolver)); 362661095bbSAlp Dener for (i = 0; i < tao->numbermonitors; i++) { 3639566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)tao->monitorcontext[i])); 36410978b7dSBarry Smith PetscCall(TaoMonitorSet(auglag->subsolver, tao->monitor[i], tao->monitorcontext[i], tao->monitordestroy[i])); 36510978b7dSBarry Smith if (tao->monitor[i] == TaoMonitorDefault || tao->monitor[i] == TaoMonitorConstraintNorm || tao->monitor[i] == TaoMonitorGlobalization || tao->monitor[i] == TaoMonitorDefaultShort) auglag->info = PETSC_TRUE; 366661095bbSAlp Dener } 3673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 368661095bbSAlp Dener } 369661095bbSAlp Dener 370661095bbSAlp Dener /* -------------------------------------------------------- */ 371661095bbSAlp Dener 372661095bbSAlp Dener /*MC 373c0298470SBarry Smith TAOALMM - Augmented Lagrangian multiplier method for solving nonlinear optimization problems with general constraints. 374661095bbSAlp Dener 375661095bbSAlp Dener Options Database Keys: 376661095bbSAlp Dener + -tao_almm_mu_init <real> - initial penalty parameter (default: 10.) 377661095bbSAlp Dener . -tao_almm_mu_factor <real> - increase factor for the penalty parameter (default: 100.) 378661095bbSAlp Dener . -tao_almm_mu_max <real> - maximum safeguard for penalty parameter updates (default: 1.e20) 379661095bbSAlp Dener . -tao_almm_mu_power_good <real> - exponential for penalty parameter when multiplier update is accepted (default: 0.9) 380661095bbSAlp Dener . -tao_almm_mu_power_bad <real> - exponential for penalty parameter when multiplier update is rejected (default: 0.1) 381661095bbSAlp Dener . -tao_almm_ye_min <real> - minimum safeguard for equality multiplier updates (default: -1.e20) 382661095bbSAlp Dener . -tao_almm_ye_max <real> - maximum safeguard for equality multiplier updates (default: 1.e20) 383661095bbSAlp Dener . -tao_almm_yi_min <real> - minimum safeguard for inequality multiplier updates (default: -1.e20) 384661095bbSAlp Dener . -tao_almm_yi_max <real> - maximum safeguard for inequality multiplier updates (default: 1.e20) 3859a09327dSAlp Dener - -tao_almm_type <phr,classic> - change formulation of the augmented Lagrangian merit function for the subproblem (default: phr) 386661095bbSAlp Dener 387661095bbSAlp Dener Level: beginner 388661095bbSAlp Dener 389661095bbSAlp Dener Notes: 390661095bbSAlp Dener This method converts a constrained problem into a sequence of unconstrained problems via the augmented 391661095bbSAlp Dener Lagrangian merit function. Bound constraints are pushed down to the subproblem without any modifications. 392661095bbSAlp Dener 393661095bbSAlp Dener Two formulations are offered for the subproblem: canonical Hestenes-Powell augmented Lagrangian with slack 394661095bbSAlp Dener variables for inequality constraints, and a slack-less Powell-Hestenes-Rockafellar (PHR) formulation utilizing a 3959a09327dSAlp Dener pointwise max() penalty on inequality constraints. The canonical augmented Lagrangian formulation may converge 3969a09327dSAlp Dener faster for smaller problems but is highly susceptible to poor step lengths in the subproblem due to the positivity 3979a09327dSAlp Dener constraint on slack variables. PHR avoids this issue by eliminating the slack variables entirely, and is highly 3989a09327dSAlp Dener desirable for problems with a large number of inequality constraints. 399661095bbSAlp Dener 4009a09327dSAlp Dener The subproblem is solved using a nested first-order TAO solver (default: `TAOBQNLS`). The user can retrieve a 4019a09327dSAlp Dener pointer to the subsolver via `TaoALMMGetSubsolver()` or pass command line arguments to it using the 4029a09327dSAlp Dener "-tao_almm_subsolver_" prefix. Currently, `TAOALMM` does not support second-order methods for the subproblem. 403661095bbSAlp Dener 404661095bbSAlp Dener .vb 405661095bbSAlp Dener while unconverged 406661095bbSAlp Dener solve argmin_x L(x) s.t. l <= x <= u 407661095bbSAlp Dener if ||c|| <= y_tol 408661095bbSAlp Dener if ||c|| <= c_tol && ||Lgrad|| <= g_tol: 409661095bbSAlp Dener problem converged, return solution 410661095bbSAlp Dener else 411661095bbSAlp Dener constraints sufficiently improved 412661095bbSAlp Dener update multipliers and tighten tolerances 413661095bbSAlp Dener endif 414661095bbSAlp Dener else 415661095bbSAlp Dener constraints did not improve 416661095bbSAlp Dener update penalty and loosen tolerances 417661095bbSAlp Dener endif 418661095bbSAlp Dener endwhile 419661095bbSAlp Dener .ve 420661095bbSAlp Dener 421c0298470SBarry Smith .seealso: `TAOALMM`, `Tao`, `TaoALMMGetType()`, `TaoALMMSetType()`, `TaoALMMSetSubsolver()`, `TaoALMMGetSubsolver()`, 422db781477SPatrick Sanan `TaoALMMGetMultipliers()`, `TaoALMMSetMultipliers()`, `TaoALMMGetPrimalIS()`, `TaoALMMGetDualIS()` 423661095bbSAlp Dener M*/ 424d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TaoCreate_ALMM(Tao tao) 425d71ae5a4SJacob Faibussowitsch { 426661095bbSAlp Dener TAO_ALMM *auglag; 427661095bbSAlp Dener 428661095bbSAlp Dener PetscFunctionBegin; 4294dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&auglag)); 430661095bbSAlp Dener 431661095bbSAlp Dener tao->ops->destroy = TaoDestroy_ALMM; 432661095bbSAlp Dener tao->ops->setup = TaoSetUp_ALMM; 433661095bbSAlp Dener tao->ops->setfromoptions = TaoSetFromOptions_ALMM; 434661095bbSAlp Dener tao->ops->view = TaoView_ALMM; 435661095bbSAlp Dener tao->ops->solve = TaoSolve_ALMM; 436661095bbSAlp Dener 437606f75f6SBarry Smith PetscCall(TaoParametersInitialize(tao)); 438606f75f6SBarry Smith PetscObjectParameterSetDefault(tao, gatol, 1.e-5); 439606f75f6SBarry Smith PetscObjectParameterSetDefault(tao, grtol, 0.0); 440606f75f6SBarry Smith PetscObjectParameterSetDefault(tao, gttol, 0.0); 441606f75f6SBarry Smith PetscObjectParameterSetDefault(tao, catol, 1.e-5); 442606f75f6SBarry Smith PetscObjectParameterSetDefault(tao, crtol, 0.0); 443661095bbSAlp Dener 444661095bbSAlp Dener tao->data = (void *)auglag; 445661095bbSAlp Dener auglag->parent = tao; 446661095bbSAlp Dener auglag->mu0 = 10.0; 447661095bbSAlp Dener auglag->mu = auglag->mu0; 448661095bbSAlp Dener auglag->mu_fac = 10.0; 449661095bbSAlp Dener auglag->mu_max = PETSC_INFINITY; 450661095bbSAlp Dener auglag->mu_pow_good = 0.9; 451661095bbSAlp Dener auglag->mu_pow_bad = 0.1; 452661095bbSAlp Dener auglag->ye_min = PETSC_NINFINITY; 453661095bbSAlp Dener auglag->ye_max = PETSC_INFINITY; 454661095bbSAlp Dener auglag->yi_min = PETSC_NINFINITY; 455661095bbSAlp Dener auglag->yi_max = PETSC_INFINITY; 456661095bbSAlp Dener auglag->ytol0 = 0.1 / PetscPowReal(auglag->mu0, auglag->mu_pow_bad); 457661095bbSAlp Dener auglag->ytol = auglag->ytol0; 458661095bbSAlp Dener auglag->gtol0 = 1.0 / auglag->mu0; 459661095bbSAlp Dener auglag->gtol = auglag->gtol0; 460661095bbSAlp Dener 461661095bbSAlp Dener auglag->sub_obj = TaoALMMComputeAugLagAndGradient_Private; 4629a09327dSAlp Dener auglag->type = TAO_ALMM_PHR; 463661095bbSAlp Dener auglag->info = PETSC_FALSE; 464661095bbSAlp Dener 4659566063dSJacob Faibussowitsch PetscCall(TaoCreate(PetscObjectComm((PetscObject)tao), &auglag->subsolver)); 4669a09327dSAlp Dener PetscCall(TaoSetType(auglag->subsolver, TAOBQNLS)); 4679566063dSJacob Faibussowitsch PetscCall(TaoSetTolerances(auglag->subsolver, auglag->gtol, 0.0, 0.0)); 4689566063dSJacob Faibussowitsch PetscCall(TaoSetMaximumIterations(auglag->subsolver, 1000)); 4699566063dSJacob Faibussowitsch PetscCall(TaoSetMaximumFunctionEvaluations(auglag->subsolver, 10000)); 4709566063dSJacob Faibussowitsch PetscCall(TaoSetFunctionLowerBound(auglag->subsolver, PETSC_NINFINITY)); 4719566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)auglag->subsolver, (PetscObject)tao, 1)); 472661095bbSAlp Dener 4739566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetType_C", TaoALMMGetType_Private)); 4749566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMSetType_C", TaoALMMSetType_Private)); 4759566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetSubsolver_C", TaoALMMGetSubsolver_Private)); 4769566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMSetSubsolver_C", TaoALMMSetSubsolver_Private)); 4779566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetMultipliers_C", TaoALMMGetMultipliers_Private)); 4789566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMSetMultipliers_C", TaoALMMSetMultipliers_Private)); 4799566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetPrimalIS_C", TaoALMMGetPrimalIS_Private)); 4809566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetDualIS_C", TaoALMMGetDualIS_Private)); 4813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 482661095bbSAlp Dener } 483661095bbSAlp Dener 484d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoALMMCombinePrimal_Private(Tao tao, Vec X, Vec S, Vec P) 485d71ae5a4SJacob Faibussowitsch { 486661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM *)tao->data; 487661095bbSAlp Dener 488661095bbSAlp Dener PetscFunctionBegin; 489661095bbSAlp Dener if (tao->ineq_constrained) { 4909566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(auglag->Pscatter[0], X, P, INSERT_VALUES, SCATTER_REVERSE)); 4919566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(auglag->Pscatter[0], X, P, INSERT_VALUES, SCATTER_REVERSE)); 4929566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(auglag->Pscatter[1], S, P, INSERT_VALUES, SCATTER_REVERSE)); 4939566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(auglag->Pscatter[1], S, P, INSERT_VALUES, SCATTER_REVERSE)); 494661095bbSAlp Dener } else { 4959566063dSJacob Faibussowitsch PetscCall(VecCopy(X, P)); 496661095bbSAlp Dener } 4973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 498661095bbSAlp Dener } 499661095bbSAlp Dener 500d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoALMMCombineDual_Private(Tao tao, Vec EQ, Vec IN, Vec Y) 501d71ae5a4SJacob Faibussowitsch { 502661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM *)tao->data; 503661095bbSAlp Dener 504661095bbSAlp Dener PetscFunctionBegin; 505661095bbSAlp Dener if (tao->eq_constrained) { 506661095bbSAlp Dener if (tao->ineq_constrained) { 5079566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(auglag->Yscatter[0], EQ, Y, INSERT_VALUES, SCATTER_REVERSE)); 5089566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(auglag->Yscatter[0], EQ, Y, INSERT_VALUES, SCATTER_REVERSE)); 5099566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(auglag->Yscatter[1], IN, Y, INSERT_VALUES, SCATTER_REVERSE)); 5109566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(auglag->Yscatter[1], IN, Y, INSERT_VALUES, SCATTER_REVERSE)); 511661095bbSAlp Dener } else { 5129566063dSJacob Faibussowitsch PetscCall(VecCopy(EQ, Y)); 513661095bbSAlp Dener } 514661095bbSAlp Dener } else { 5159566063dSJacob Faibussowitsch PetscCall(VecCopy(IN, Y)); 516661095bbSAlp Dener } 5173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 518661095bbSAlp Dener } 519661095bbSAlp Dener 520d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoALMMSplitPrimal_Private(Tao tao, Vec P, Vec X, Vec S) 521d71ae5a4SJacob Faibussowitsch { 522661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM *)tao->data; 523661095bbSAlp Dener 524661095bbSAlp Dener PetscFunctionBegin; 525661095bbSAlp Dener if (tao->ineq_constrained) { 5269566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(auglag->Pscatter[0], P, X, INSERT_VALUES, SCATTER_FORWARD)); 5279566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(auglag->Pscatter[0], P, X, INSERT_VALUES, SCATTER_FORWARD)); 5289566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(auglag->Pscatter[1], P, S, INSERT_VALUES, SCATTER_FORWARD)); 5299566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(auglag->Pscatter[1], P, S, INSERT_VALUES, SCATTER_FORWARD)); 530661095bbSAlp Dener } else { 5319566063dSJacob Faibussowitsch PetscCall(VecCopy(P, X)); 532661095bbSAlp Dener } 5333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 534661095bbSAlp Dener } 535661095bbSAlp Dener 536661095bbSAlp Dener /* this assumes that the latest constraints are stored in Ce and Ci, and also combined in C */ 537d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoALMMComputeOptimalityNorms_Private(Tao tao) 538d71ae5a4SJacob Faibussowitsch { 539661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM *)tao->data; 540661095bbSAlp Dener 541661095bbSAlp Dener PetscFunctionBegin; 542661095bbSAlp Dener /* if bounded, project the gradient */ 5431baa6e33SBarry Smith if (tao->bounded) PetscCall(VecBoundGradientProjection(auglag->LgradX, auglag->Px, tao->XL, tao->XU, auglag->LgradX)); 544661095bbSAlp Dener if (auglag->type == TAO_ALMM_PHR) { 5459566063dSJacob Faibussowitsch PetscCall(VecNorm(auglag->LgradX, NORM_INFINITY, &auglag->gnorm)); 546661095bbSAlp Dener auglag->cenorm = 0.0; 54748a46eb9SPierre Jolivet if (tao->eq_constrained) PetscCall(VecNorm(auglag->Ce, NORM_INFINITY, &auglag->cenorm)); 548661095bbSAlp Dener auglag->cinorm = 0.0; 549661095bbSAlp Dener if (tao->ineq_constrained) { 5509566063dSJacob Faibussowitsch PetscCall(VecCopy(auglag->Yi, auglag->Ciwork)); 5519566063dSJacob Faibussowitsch PetscCall(VecScale(auglag->Ciwork, -1.0 / auglag->mu)); 5529566063dSJacob Faibussowitsch PetscCall(VecPointwiseMax(auglag->Ciwork, auglag->Ci, auglag->Ciwork)); 5539566063dSJacob Faibussowitsch PetscCall(VecNorm(auglag->Ciwork, NORM_INFINITY, &auglag->cinorm)); 554661095bbSAlp Dener } 555661095bbSAlp Dener auglag->cnorm_old = auglag->cnorm; 556661095bbSAlp Dener auglag->cnorm = PetscMax(auglag->cenorm, auglag->cinorm); 557661095bbSAlp Dener auglag->ytol = auglag->ytol0 * auglag->cnorm_old; 558661095bbSAlp Dener } else { 5599566063dSJacob Faibussowitsch PetscCall(VecNorm(auglag->LgradX, NORM_2, &auglag->gnorm)); 5609566063dSJacob Faibussowitsch PetscCall(VecNorm(auglag->C, NORM_2, &auglag->cnorm)); 561661095bbSAlp Dener } 5623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 563661095bbSAlp Dener } 564661095bbSAlp Dener 5659a09327dSAlp Dener static PetscErrorCode TaoALMMEvaluateIterate_Private(Tao tao) 566d71ae5a4SJacob Faibussowitsch { 567661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM *)tao->data; 568661095bbSAlp Dener 569661095bbSAlp Dener PetscFunctionBegin; 570661095bbSAlp Dener /* split solution into primal and slack components */ 5719566063dSJacob Faibussowitsch PetscCall(TaoALMMSplitPrimal_Private(tao, auglag->P, auglag->Px, auglag->Ps)); 572661095bbSAlp Dener 573661095bbSAlp Dener /* compute f, df/dx and the constraints */ 5749566063dSJacob Faibussowitsch PetscCall(TaoComputeObjectiveAndGradient(tao, auglag->Px, &auglag->fval, auglag->LgradX)); 575661095bbSAlp Dener if (tao->eq_constrained) { 5769566063dSJacob Faibussowitsch PetscCall(TaoComputeEqualityConstraints(tao, auglag->Px, auglag->Ce)); 5779566063dSJacob Faibussowitsch PetscCall(TaoComputeJacobianEquality(tao, auglag->Px, auglag->Ae, auglag->Ae)); 578661095bbSAlp Dener } 579661095bbSAlp Dener if (tao->ineq_constrained) { 5809566063dSJacob Faibussowitsch PetscCall(TaoComputeInequalityConstraints(tao, auglag->Px, auglag->Ci)); 5819566063dSJacob Faibussowitsch PetscCall(TaoComputeJacobianInequality(tao, auglag->Px, auglag->Ai, auglag->Ai)); 582661095bbSAlp Dener switch (auglag->type) { 5837721a36fSStefano Zampini case TAO_ALMM_CLASSIC: 584661095bbSAlp Dener /* classic formulation converts inequality to equality constraints via slack variables */ 5859566063dSJacob Faibussowitsch PetscCall(VecAXPY(auglag->Ci, -1.0, auglag->Ps)); 586661095bbSAlp Dener break; 5877721a36fSStefano Zampini case TAO_ALMM_PHR: 588661095bbSAlp Dener /* PHR is based on Ci <= 0 while TAO defines Ci >= 0 so we hit it with a negative sign */ 5899566063dSJacob Faibussowitsch PetscCall(VecScale(auglag->Ci, -1.0)); 5909566063dSJacob Faibussowitsch PetscCall(MatScale(auglag->Ai, -1.0)); 591661095bbSAlp Dener break; 592d71ae5a4SJacob Faibussowitsch default: 593d71ae5a4SJacob Faibussowitsch break; 594661095bbSAlp Dener } 595661095bbSAlp Dener } 596661095bbSAlp Dener /* combine constraints into one vector */ 5979566063dSJacob Faibussowitsch PetscCall(TaoALMMCombineDual_Private(tao, auglag->Ce, auglag->Ci, auglag->C)); 5983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 599661095bbSAlp Dener } 600661095bbSAlp Dener 601661095bbSAlp Dener /* 602661095bbSAlp Dener Lphr = f + 0.5*mu*[ (Ce + Ye/mu)^T (Ce + Ye/mu) + pmin(0, Ci + Yi/mu)^T pmin(0, Ci + Yi/mu)] 603661095bbSAlp Dener 604661095bbSAlp Dener dLphr/dX = dF/dX + mu*[ (Ce + Ye/mu)^T Ae + pmin(0, Ci + Yi/mu)^T Ai] 605661095bbSAlp Dener 606661095bbSAlp Dener dLphr/dS = 0 607661095bbSAlp Dener */ 608d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoALMMComputePHRLagAndGradient_Private(Tao tao) 609d71ae5a4SJacob Faibussowitsch { 610661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM *)tao->data; 611661095bbSAlp Dener PetscReal eq_norm = 0.0, ineq_norm = 0.0; 612661095bbSAlp Dener 613661095bbSAlp Dener PetscFunctionBegin; 6149a09327dSAlp Dener PetscCall(TaoALMMEvaluateIterate_Private(tao)); 615661095bbSAlp Dener if (tao->eq_constrained) { 616661095bbSAlp Dener /* Ce_work = mu*(Ce + Ye/mu) */ 6179566063dSJacob Faibussowitsch PetscCall(VecWAXPY(auglag->Cework, 1.0 / auglag->mu, auglag->Ye, auglag->Ce)); 6189566063dSJacob Faibussowitsch PetscCall(VecDot(auglag->Cework, auglag->Cework, &eq_norm)); /* contribution to scalar Lagrangian */ 6199566063dSJacob Faibussowitsch PetscCall(VecScale(auglag->Cework, auglag->mu)); 620661095bbSAlp Dener /* dL/dX += mu*(Ce + Ye/mu)^T Ae */ 6219566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(auglag->Ae, auglag->Cework, auglag->LgradX, auglag->LgradX)); 622661095bbSAlp Dener } 623661095bbSAlp Dener if (tao->ineq_constrained) { 624661095bbSAlp Dener /* Ci_work = mu * pmax(0, Ci + Yi/mu) where pmax() is pointwise max() */ 6259566063dSJacob Faibussowitsch PetscCall(VecWAXPY(auglag->Ciwork, 1.0 / auglag->mu, auglag->Yi, auglag->Ci)); 6269566063dSJacob Faibussowitsch PetscCall(VecPointwiseMax(auglag->Ciwork, auglag->Cizero, auglag->Ciwork)); 6279566063dSJacob Faibussowitsch PetscCall(VecDot(auglag->Ciwork, auglag->Ciwork, &ineq_norm)); /* contribution to scalar Lagrangian */ 628661095bbSAlp Dener /* dL/dX += mu * pmax(0, Ci + Yi/mu)^T Ai */ 6299566063dSJacob Faibussowitsch PetscCall(VecScale(auglag->Ciwork, auglag->mu)); 6309566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(auglag->Ai, auglag->Ciwork, auglag->LgradX, auglag->LgradX)); 631a5b23f4aSJose E. Roman /* dL/dS = 0 because there are no slacks in PHR */ 6329566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(auglag->LgradS)); 633661095bbSAlp Dener } 634661095bbSAlp Dener /* combine gradient together */ 6359566063dSJacob Faibussowitsch PetscCall(TaoALMMCombinePrimal_Private(tao, auglag->LgradX, auglag->LgradS, auglag->G)); 636661095bbSAlp 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)] */ 6375f80ce2aSJacob Faibussowitsch auglag->Lval = auglag->fval + 0.5 * auglag->mu * (eq_norm + ineq_norm); 6383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 639661095bbSAlp Dener } 640661095bbSAlp Dener 641661095bbSAlp Dener /* 642661095bbSAlp Dener Lc = F + Ye^TCe + Yi^T(Ci - S) + 0.5*mu*[Ce^TCe + (Ci - S)^T(Ci - S)] 643661095bbSAlp Dener 644661095bbSAlp Dener dLc/dX = dF/dX + Ye^TAe + Yi^TAi + 0.5*mu*[Ce^TAe + (Ci - S)^TAi] 645661095bbSAlp Dener 646661095bbSAlp Dener dLc/dS = -[Yi + mu*(Ci - S)] 647661095bbSAlp Dener */ 648d71ae5a4SJacob Faibussowitsch static PetscErrorCode TaoALMMComputeAugLagAndGradient_Private(Tao tao) 649d71ae5a4SJacob Faibussowitsch { 650661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM *)tao->data; 651661095bbSAlp Dener PetscReal yeTce = 0.0, yiTcims = 0.0, ceTce = 0.0, cimsTcims = 0.0; 652661095bbSAlp Dener 653661095bbSAlp Dener PetscFunctionBegin; 6549a09327dSAlp Dener PetscCall(TaoALMMEvaluateIterate_Private(tao)); 655661095bbSAlp Dener if (tao->eq_constrained) { 656661095bbSAlp Dener /* compute scalar contributions */ 6579566063dSJacob Faibussowitsch PetscCall(VecDot(auglag->Ye, auglag->Ce, &yeTce)); 6589566063dSJacob Faibussowitsch PetscCall(VecDot(auglag->Ce, auglag->Ce, &ceTce)); 659661095bbSAlp Dener /* dL/dX += ye^T Ae */ 6609566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(auglag->Ae, auglag->Ye, auglag->LgradX, auglag->LgradX)); 6619a09327dSAlp Dener /* dL/dX += 0.5 * mu * ce^T Ae */ 6629566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(auglag->Ae, auglag->Ce, auglag->Xwork)); 6639a09327dSAlp Dener PetscCall(VecAXPY(auglag->LgradX, 0.5 * auglag->mu, auglag->Xwork)); 664661095bbSAlp Dener } 665661095bbSAlp Dener if (tao->ineq_constrained) { 666661095bbSAlp Dener /* compute scalar contributions */ 6679566063dSJacob Faibussowitsch PetscCall(VecDot(auglag->Yi, auglag->Ci, &yiTcims)); 6689566063dSJacob Faibussowitsch PetscCall(VecDot(auglag->Ci, auglag->Ci, &cimsTcims)); 669661095bbSAlp Dener /* dL/dX += yi^T Ai */ 6709566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(auglag->Ai, auglag->Yi, auglag->LgradX, auglag->LgradX)); 6719a09327dSAlp Dener /* dL/dX += 0.5 * mu * (ci - s)^T Ai */ 6729566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(auglag->Ai, auglag->Ci, auglag->Xwork)); 6739a09327dSAlp Dener PetscCall(VecAXPY(auglag->LgradX, 0.5 * auglag->mu, auglag->Xwork)); 674661095bbSAlp Dener /* dL/dS = -[yi + mu*(ci - s)] */ 6759566063dSJacob Faibussowitsch PetscCall(VecWAXPY(auglag->LgradS, auglag->mu, auglag->Ci, auglag->Yi)); 6769566063dSJacob Faibussowitsch PetscCall(VecScale(auglag->LgradS, -1.0)); 677661095bbSAlp Dener } 678661095bbSAlp Dener /* combine gradient together */ 6799566063dSJacob Faibussowitsch PetscCall(TaoALMMCombinePrimal_Private(tao, auglag->LgradX, auglag->LgradS, auglag->G)); 680661095bbSAlp Dener /* compute L = f + ye^T ce + yi^T (ci - s) + 0.5*mu*||ce||^2 + 0.5*mu*||ci - s||^2 */ 681661095bbSAlp Dener auglag->Lval = auglag->fval + yeTce + yiTcims + 0.5 * auglag->mu * (ceTce + cimsTcims); 6823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 683661095bbSAlp Dener } 684661095bbSAlp Dener 685d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoALMMSubsolverObjective_Private(Tao tao, Vec P, PetscReal *Lval, void *ctx) 686d71ae5a4SJacob Faibussowitsch { 6877721a36fSStefano Zampini TAO_ALMM *auglag = (TAO_ALMM *)ctx; 6887721a36fSStefano Zampini 6897721a36fSStefano Zampini PetscFunctionBegin; 6909566063dSJacob Faibussowitsch PetscCall(VecCopy(P, auglag->P)); 6919566063dSJacob Faibussowitsch PetscCall((*auglag->sub_obj)(auglag->parent)); 6927721a36fSStefano Zampini *Lval = auglag->Lval; 6933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6947721a36fSStefano Zampini } 6957721a36fSStefano Zampini 696d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoALMMSubsolverObjectiveAndGradient_Private(Tao tao, Vec P, PetscReal *Lval, Vec G, void *ctx) 697d71ae5a4SJacob Faibussowitsch { 698661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM *)ctx; 699661095bbSAlp Dener 700661095bbSAlp Dener PetscFunctionBegin; 7019566063dSJacob Faibussowitsch PetscCall(VecCopy(P, auglag->P)); 7029566063dSJacob Faibussowitsch PetscCall((*auglag->sub_obj)(auglag->parent)); 7039566063dSJacob Faibussowitsch PetscCall(VecCopy(auglag->G, G)); 704661095bbSAlp Dener *Lval = auglag->Lval; 7053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 706661095bbSAlp Dener } 707