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 13661095bbSAlp Dener static PetscErrorCode TaoSolve_ALMM(Tao tao) 14661095bbSAlp Dener { 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)); 37*dbbe0bcdSBarry Smith PetscUseTypeMethod(tao,convergencetest , tao->cnvP); 38661095bbSAlp Dener /* set initial penalty factor and inner solver tolerance */ 39661095bbSAlp Dener switch (auglag->type) { 407721a36fSStefano Zampini case TAO_ALMM_CLASSIC: 41661095bbSAlp Dener auglag->mu = auglag->mu0; 42661095bbSAlp Dener break; 437721a36fSStefano Zampini case TAO_ALMM_PHR: 44661095bbSAlp Dener auglag->cenorm = 0.0; 45661095bbSAlp Dener if (tao->eq_constrained) { 469566063dSJacob Faibussowitsch PetscCall(VecDot(auglag->Ce, auglag->Ce, &auglag->cenorm)); 47661095bbSAlp Dener } 48661095bbSAlp Dener auglag->cinorm = 0.0; 49661095bbSAlp Dener if (tao->ineq_constrained) { 509566063dSJacob Faibussowitsch PetscCall(VecCopy(auglag->Ci, auglag->Ciwork)); 519566063dSJacob Faibussowitsch PetscCall(VecScale(auglag->Ciwork, -1.0)); 529566063dSJacob Faibussowitsch PetscCall(VecPointwiseMax(auglag->Ciwork, auglag->Cizero, auglag->Ciwork)); 539566063dSJacob Faibussowitsch PetscCall(VecDot(auglag->Ciwork, auglag->Ciwork, &auglag->cinorm)); 54661095bbSAlp Dener } 55661095bbSAlp Dener /* determine initial penalty factor based on the balance of constraint violation and objective function value */ 56661095bbSAlp Dener auglag->mu = PetscMax(1.e-6, PetscMin(10.0, 2.0*PetscAbsReal(auglag->fval)/(auglag->cenorm + auglag->cinorm))); 57661095bbSAlp Dener break; 58661095bbSAlp Dener default: 59661095bbSAlp Dener break; 60661095bbSAlp Dener } 61661095bbSAlp Dener auglag->gtol = auglag->gtol0; 6263a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(tao,"Initial penalty: %.2g\n",(double)auglag->mu)); 63661095bbSAlp Dener 64661095bbSAlp Dener /* start aug-lag outer loop */ 65661095bbSAlp Dener while (tao->reason == TAO_CONTINUE_ITERATING) { 66661095bbSAlp Dener ++tao->niter; 67661095bbSAlp Dener /* update subsolver tolerance */ 6863a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(tao,"Subsolver tolerance: ||G|| <= %e\n",(double)auglag->gtol)); 699566063dSJacob Faibussowitsch PetscCall(TaoSetTolerances(auglag->subsolver, auglag->gtol, 0.0, 0.0)); 70661095bbSAlp Dener /* solve the bound-constrained or unconstrained subproblem */ 719566063dSJacob Faibussowitsch PetscCall(TaoSolve(auglag->subsolver)); 729566063dSJacob Faibussowitsch PetscCall(TaoGetConvergedReason(auglag->subsolver, &reason)); 73661095bbSAlp Dener tao->ksp_its += auglag->subsolver->ksp_its; 74661095bbSAlp Dener if (reason != TAO_CONVERGED_GATOL) { 759566063dSJacob Faibussowitsch PetscCall(PetscInfo(tao,"Subsolver failed to converge, reason: %s\n",TaoConvergedReasons[reason])); 76661095bbSAlp Dener } 77661095bbSAlp Dener /* evaluate solution and test convergence */ 789566063dSJacob Faibussowitsch PetscCall((*auglag->sub_obj)(tao)); 799566063dSJacob Faibussowitsch PetscCall(TaoALMMComputeOptimalityNorms_Private(tao)); 80661095bbSAlp Dener /* decide whether to update multipliers or not */ 81661095bbSAlp Dener updated = 0.0; 82661095bbSAlp Dener if (auglag->cnorm <= auglag->ytol) { 8363a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(tao,"Multipliers updated: ||C|| <= %e\n",(double)auglag->ytol)); 84661095bbSAlp Dener /* constraints are good, update multipliers and convergence tolerances */ 85661095bbSAlp Dener if (tao->eq_constrained) { 869566063dSJacob Faibussowitsch PetscCall(VecAXPY(auglag->Ye, auglag->mu, auglag->Ce)); 879566063dSJacob Faibussowitsch PetscCall(VecSet(auglag->Cework, auglag->ye_max)); 889566063dSJacob Faibussowitsch PetscCall(VecPointwiseMin(auglag->Ye, auglag->Cework, auglag->Ye)); 899566063dSJacob Faibussowitsch PetscCall(VecSet(auglag->Cework, auglag->ye_min)); 909566063dSJacob Faibussowitsch PetscCall(VecPointwiseMax(auglag->Ye, auglag->Cework, auglag->Ye)); 91661095bbSAlp Dener } 92661095bbSAlp Dener if (tao->ineq_constrained) { 939566063dSJacob Faibussowitsch PetscCall(VecAXPY(auglag->Yi, auglag->mu, auglag->Ci)); 949566063dSJacob Faibussowitsch PetscCall(VecSet(auglag->Ciwork, auglag->yi_max)); 959566063dSJacob Faibussowitsch PetscCall(VecPointwiseMin(auglag->Yi, auglag->Ciwork, auglag->Yi)); 969566063dSJacob Faibussowitsch PetscCall(VecSet(auglag->Ciwork, auglag->yi_min)); 979566063dSJacob Faibussowitsch PetscCall(VecPointwiseMax(auglag->Yi, auglag->Ciwork, auglag->Yi)); 98661095bbSAlp Dener } 99661095bbSAlp Dener /* tolerances are updated only for non-PHR methods */ 100661095bbSAlp Dener if (auglag->type != TAO_ALMM_PHR) { 101661095bbSAlp Dener auglag->ytol = PetscMax(tao->catol, auglag->ytol/PetscPowReal(auglag->mu, auglag->mu_pow_good)); 102661095bbSAlp Dener auglag->gtol = PetscMax(tao->gatol, auglag->gtol/auglag->mu); 103661095bbSAlp Dener } 104661095bbSAlp Dener updated = 1.0; 105661095bbSAlp Dener } else { 106661095bbSAlp Dener /* constraints are bad, update penalty factor */ 107661095bbSAlp Dener auglag->mu = PetscMin(auglag->mu_max, auglag->mu_fac*auglag->mu); 108661095bbSAlp Dener /* tolerances are reset only for non-PHR methods */ 109661095bbSAlp Dener if (auglag->type != TAO_ALMM_PHR) { 110661095bbSAlp Dener auglag->ytol = PetscMax(tao->catol, 0.1/PetscPowReal(auglag->mu, auglag->mu_pow_bad)); 111661095bbSAlp Dener auglag->gtol = PetscMax(tao->gatol, 1.0/auglag->mu); 112661095bbSAlp Dener } 11363a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(tao,"Penalty increased: mu = %.2g\n",(double)auglag->mu)); 114661095bbSAlp Dener } 1159566063dSJacob Faibussowitsch PetscCall(TaoLogConvergenceHistory(tao, auglag->fval, auglag->gnorm, auglag->cnorm, tao->ksp_its)); 1169566063dSJacob Faibussowitsch PetscCall(TaoMonitor(tao, tao->niter, auglag->fval, auglag->gnorm, auglag->cnorm, updated)); 117*dbbe0bcdSBarry Smith PetscUseTypeMethod(tao,convergencetest , tao->cnvP); 118661095bbSAlp Dener } 119661095bbSAlp Dener 120661095bbSAlp Dener PetscFunctionReturn(0); 121661095bbSAlp Dener } 122661095bbSAlp Dener 123661095bbSAlp Dener static PetscErrorCode TaoView_ALMM(Tao tao,PetscViewer viewer) 124661095bbSAlp Dener { 125661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM*)tao->data; 126661095bbSAlp Dener PetscBool isascii; 127661095bbSAlp Dener 128661095bbSAlp Dener PetscFunctionBegin; 1299566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 1309566063dSJacob Faibussowitsch PetscCall(TaoView(auglag->subsolver,viewer)); 1319566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 1329566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 133661095bbSAlp Dener if (isascii) { 1349566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 1359566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "ALMM Formulation Type: %s\n", TaoALMMTypes[auglag->type])); 1369566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 137661095bbSAlp Dener } 138661095bbSAlp Dener PetscFunctionReturn(0); 139661095bbSAlp Dener } 140661095bbSAlp Dener 141661095bbSAlp Dener static PetscErrorCode TaoSetUp_ALMM(Tao tao) 142661095bbSAlp Dener { 143661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM*)tao->data; 144661095bbSAlp Dener VecType vec_type; 145661095bbSAlp Dener Vec SL, SU; 146661095bbSAlp Dener PetscBool is_cg = PETSC_FALSE, is_lmvm = PETSC_FALSE; 147661095bbSAlp Dener 148661095bbSAlp Dener PetscFunctionBegin; 1493c859ba3SBarry 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."); 1503c859ba3SBarry 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."); 1519566063dSJacob Faibussowitsch PetscCall(TaoComputeVariableBounds(tao)); 152661095bbSAlp Dener /* alias base vectors and create extras */ 1539566063dSJacob Faibussowitsch PetscCall(VecGetType(tao->solution, &vec_type)); 154661095bbSAlp Dener auglag->Px = tao->solution; 155661095bbSAlp Dener if (!tao->gradient) { /* base gradient */ 1569566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &tao->gradient)); 157661095bbSAlp Dener } 158661095bbSAlp Dener auglag->LgradX = tao->gradient; 159661095bbSAlp Dener if (!auglag->Xwork) { /* opt var work vector */ 1609566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tao->solution, &auglag->Xwork)); 161661095bbSAlp Dener } 162661095bbSAlp Dener if (tao->eq_constrained) { 163661095bbSAlp Dener auglag->Ce = tao->constraints_equality; 164661095bbSAlp Dener auglag->Ae = tao->jacobian_equality; 165661095bbSAlp Dener if (!auglag->Ye) { /* equality multipliers */ 1669566063dSJacob Faibussowitsch PetscCall(VecDuplicate(auglag->Ce, &auglag->Ye)); 167661095bbSAlp Dener } 168661095bbSAlp Dener if (!auglag->Cework) { 1699566063dSJacob Faibussowitsch PetscCall(VecDuplicate(auglag->Ce, &auglag->Cework)); 170661095bbSAlp Dener } 171661095bbSAlp Dener } 172661095bbSAlp Dener if (tao->ineq_constrained) { 173661095bbSAlp Dener auglag->Ci = tao->constraints_inequality; 174661095bbSAlp Dener auglag->Ai = tao->jacobian_inequality; 175661095bbSAlp Dener if (!auglag->Yi) { /* inequality multipliers */ 1769566063dSJacob Faibussowitsch PetscCall(VecDuplicate(auglag->Ci, &auglag->Yi)); 177661095bbSAlp Dener } 178661095bbSAlp Dener if (!auglag->Ciwork) { 1799566063dSJacob Faibussowitsch PetscCall(VecDuplicate(auglag->Ci, &auglag->Ciwork)); 180661095bbSAlp Dener } 181661095bbSAlp Dener if (!auglag->Cizero) { 1829566063dSJacob Faibussowitsch PetscCall(VecDuplicate(auglag->Ci, &auglag->Cizero)); 1839566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(auglag->Cizero)); 184661095bbSAlp Dener } 185661095bbSAlp Dener if (!auglag->Ps) { /* slack vars */ 1869566063dSJacob Faibussowitsch PetscCall(VecDuplicate(auglag->Ci, &auglag->Ps)); 187661095bbSAlp Dener } 188661095bbSAlp Dener if (!auglag->LgradS) { /* slack component of Lagrangian gradient */ 1899566063dSJacob Faibussowitsch PetscCall(VecDuplicate(auglag->Ci, &auglag->LgradS)); 190661095bbSAlp Dener } 191661095bbSAlp Dener /* create vector for combined primal space and the associated communication objects */ 192661095bbSAlp Dener if (!auglag->P) { 1939566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2, &auglag->Parr)); 194661095bbSAlp Dener auglag->Parr[0] = auglag->Px; auglag->Parr[1] = auglag->Ps; 1959566063dSJacob Faibussowitsch PetscCall(VecConcatenate(2, auglag->Parr, &auglag->P, &auglag->Pis)); 1969566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2, &auglag->Pscatter)); 1979566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(auglag->P, auglag->Pis[0], auglag->Px, NULL, &auglag->Pscatter[0])); 1989566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(auglag->P, auglag->Pis[1], auglag->Ps, NULL, &auglag->Pscatter[1])); 199661095bbSAlp Dener } 200661095bbSAlp Dener if (tao->eq_constrained) { 201661095bbSAlp Dener /* create vector for combined dual space and the associated communication objects */ 202661095bbSAlp Dener if (!auglag->Y) { 2039566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2, &auglag->Yarr)); 204661095bbSAlp Dener auglag->Yarr[0] = auglag->Ye; auglag->Yarr[1] = auglag->Yi; 2059566063dSJacob Faibussowitsch PetscCall(VecConcatenate(2, auglag->Yarr, &auglag->Y, &auglag->Yis)); 2069566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2, &auglag->Yscatter)); 2079566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(auglag->Y, auglag->Yis[0], auglag->Ye, NULL, &auglag->Yscatter[0])); 2089566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(auglag->Y, auglag->Yis[1], auglag->Yi, NULL, &auglag->Yscatter[1])); 209661095bbSAlp Dener } 210661095bbSAlp Dener if (!auglag->C) { 2119566063dSJacob Faibussowitsch PetscCall(VecDuplicate(auglag->Y, &auglag->C)); 212661095bbSAlp Dener } 213661095bbSAlp Dener } else { 214661095bbSAlp Dener if (!auglag->C) { 215661095bbSAlp Dener auglag->C = auglag->Ci; 216661095bbSAlp Dener } 217661095bbSAlp Dener if (!auglag->Y) { 218661095bbSAlp Dener auglag->Y = auglag->Yi; 219661095bbSAlp Dener } 220661095bbSAlp Dener } 221661095bbSAlp Dener } else { 222661095bbSAlp Dener if (!auglag->P) { 223661095bbSAlp Dener auglag->P = auglag->Px; 224661095bbSAlp Dener } 225661095bbSAlp Dener if (!auglag->G) { 226661095bbSAlp Dener auglag->G = auglag->LgradX; 227661095bbSAlp Dener } 228661095bbSAlp Dener if (!auglag->C) { 229661095bbSAlp Dener auglag->C = auglag->Ce; 230661095bbSAlp Dener } 231661095bbSAlp Dener if (!auglag->Y) { 232661095bbSAlp Dener auglag->Y = auglag->Ye; 233661095bbSAlp Dener } 234661095bbSAlp Dener } 235661095bbSAlp Dener /* initialize parameters */ 236661095bbSAlp Dener if (auglag->type == TAO_ALMM_PHR) { 237661095bbSAlp Dener auglag->mu_fac = 10.0; 238661095bbSAlp Dener auglag->yi_min = 0.0; 239661095bbSAlp Dener auglag->ytol0 = 0.5; 240661095bbSAlp Dener auglag->gtol0 = tao->gatol; 241661095bbSAlp Dener if (tao->gatol_changed && tao->catol_changed) { 2429566063dSJacob Faibussowitsch PetscCall(PetscInfo(tao,"TAOALMM with PHR: different gradient and constraint tolerances are not supported, setting catol = gatol\n")); 243661095bbSAlp Dener tao->catol = tao->gatol; 244661095bbSAlp Dener } 245661095bbSAlp Dener } 246661095bbSAlp Dener /* set the Lagrangian formulation type for the subsolver */ 247661095bbSAlp Dener switch (auglag->type) { 2487721a36fSStefano Zampini case TAO_ALMM_CLASSIC: 249661095bbSAlp Dener auglag->sub_obj = TaoALMMComputeAugLagAndGradient_Private; 250661095bbSAlp Dener break; 2517721a36fSStefano Zampini case TAO_ALMM_PHR: 252661095bbSAlp Dener auglag->sub_obj = TaoALMMComputePHRLagAndGradient_Private; 253661095bbSAlp Dener break; 254661095bbSAlp Dener default: 255661095bbSAlp Dener break; 256661095bbSAlp Dener } 257661095bbSAlp Dener /* set up the subsolver */ 2589566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(auglag->subsolver, auglag->P)); 2599566063dSJacob Faibussowitsch PetscCall(TaoSetObjective(auglag->subsolver, TaoALMMSubsolverObjective_Private, (void*)auglag)); 2609566063dSJacob Faibussowitsch PetscCall(TaoSetObjectiveAndGradient(auglag->subsolver, NULL, TaoALMMSubsolverObjectiveAndGradient_Private, (void*)auglag)); 261661095bbSAlp Dener if (tao->bounded) { 262661095bbSAlp Dener /* make sure that the subsolver is a bound-constrained method */ 2639566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)auglag->subsolver, TAOCG, &is_cg)); 2649566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)auglag->subsolver, TAOLMVM, &is_lmvm)); 265661095bbSAlp Dener if (is_cg) { 2669566063dSJacob Faibussowitsch PetscCall(TaoSetType(auglag->subsolver, TAOBNCG)); 2679566063dSJacob Faibussowitsch PetscCall(PetscInfo(tao,"TAOCG detected for bound-constrained problem, switching to TAOBNCG instead.")); 268661095bbSAlp Dener } 269661095bbSAlp Dener if (is_lmvm) { 2709566063dSJacob Faibussowitsch PetscCall(TaoSetType(auglag->subsolver, TAOBQNLS)); 2719566063dSJacob Faibussowitsch PetscCall(PetscInfo(tao,"TAOLMVM detected for bound-constrained problem, switching to TAOBQNLS instead.")); 272661095bbSAlp Dener } 273661095bbSAlp Dener /* create lower and upper bound clone vectors for subsolver */ 274661095bbSAlp Dener if (!auglag->PL) { 2759566063dSJacob Faibussowitsch PetscCall(VecDuplicate(auglag->P, &auglag->PL)); 276661095bbSAlp Dener } 277661095bbSAlp Dener if (!auglag->PU) { 2789566063dSJacob Faibussowitsch PetscCall(VecDuplicate(auglag->P, &auglag->PU)); 279661095bbSAlp Dener } 280661095bbSAlp Dener if (tao->ineq_constrained) { 281661095bbSAlp Dener /* create lower and upper bounds for slack, set lower to 0 */ 2829566063dSJacob Faibussowitsch PetscCall(VecDuplicate(auglag->Ci, &SL)); 2839566063dSJacob Faibussowitsch PetscCall(VecSet(SL, 0.0)); 2849566063dSJacob Faibussowitsch PetscCall(VecDuplicate(auglag->Ci, &SU)); 2859566063dSJacob Faibussowitsch PetscCall(VecSet(SU, PETSC_INFINITY)); 286661095bbSAlp Dener /* combine opt var bounds with slack bounds */ 2879566063dSJacob Faibussowitsch PetscCall(TaoALMMCombinePrimal_Private(tao, tao->XL, SL, auglag->PL)); 2889566063dSJacob Faibussowitsch PetscCall(TaoALMMCombinePrimal_Private(tao, tao->XU, SU, auglag->PU)); 289661095bbSAlp Dener /* destroy work vectors */ 2909566063dSJacob Faibussowitsch PetscCall(VecDestroy(&SL)); 2919566063dSJacob Faibussowitsch PetscCall(VecDestroy(&SU)); 292661095bbSAlp Dener } else { 293661095bbSAlp Dener /* no inequality constraints, just copy bounds into the subsolver */ 2949566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->XL, auglag->PL)); 2959566063dSJacob Faibussowitsch PetscCall(VecCopy(tao->XU, auglag->PU)); 296661095bbSAlp Dener } 2979566063dSJacob Faibussowitsch PetscCall(TaoSetVariableBounds(auglag->subsolver, auglag->PL, auglag->PU)); 298661095bbSAlp Dener } 2999566063dSJacob Faibussowitsch PetscCall(TaoSetUp(auglag->subsolver)); 300661095bbSAlp Dener auglag->G = auglag->subsolver->gradient; 301661095bbSAlp Dener 302661095bbSAlp Dener PetscFunctionReturn(0); 303661095bbSAlp Dener } 304661095bbSAlp Dener 305661095bbSAlp Dener static PetscErrorCode TaoDestroy_ALMM(Tao tao) 306661095bbSAlp Dener { 307661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM*)tao->data; 308661095bbSAlp Dener 309661095bbSAlp Dener PetscFunctionBegin; 3109566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&auglag->subsolver)); 311661095bbSAlp Dener if (tao->setupcalled) { 3129566063dSJacob Faibussowitsch PetscCall(VecDestroy(&auglag->Xwork)); /* opt work */ 313661095bbSAlp Dener if (tao->eq_constrained) { 3149566063dSJacob Faibussowitsch PetscCall(VecDestroy(&auglag->Ye)); /* equality multipliers */ 3159566063dSJacob Faibussowitsch PetscCall(VecDestroy(&auglag->Cework)); /* equality work vector */ 316661095bbSAlp Dener } 317661095bbSAlp Dener if (tao->ineq_constrained) { 3189566063dSJacob Faibussowitsch PetscCall(VecDestroy(&auglag->Ps)); /* slack vars */ 319661095bbSAlp Dener auglag->Parr[0] = NULL; /* clear pointer to tao->solution, will be destroyed by TaoDestroy() shell */ 3209566063dSJacob Faibussowitsch PetscCall(PetscFree(auglag->Parr)); /* array of primal vectors */ 3219566063dSJacob Faibussowitsch PetscCall(VecDestroy(&auglag->LgradS)); /* slack grad */ 3229566063dSJacob Faibussowitsch PetscCall(VecDestroy(&auglag->Cizero)); /* zero vector for pointwise max */ 3239566063dSJacob Faibussowitsch PetscCall(VecDestroy(&auglag->Yi)); /* inequality multipliers */ 3249566063dSJacob Faibussowitsch PetscCall(VecDestroy(&auglag->Ciwork)); /* inequality work vector */ 3259566063dSJacob Faibussowitsch PetscCall(VecDestroy(&auglag->P)); /* combo primal */ 3269566063dSJacob Faibussowitsch PetscCall(ISDestroy(&auglag->Pis[0])); /* index set for X inside P */ 3279566063dSJacob Faibussowitsch PetscCall(ISDestroy(&auglag->Pis[1])); /* index set for S inside P */ 3289566063dSJacob Faibussowitsch PetscCall(PetscFree(auglag->Pis)); /* array of P index sets */ 3299566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&auglag->Pscatter[0])); 3309566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&auglag->Pscatter[1])); 3319566063dSJacob Faibussowitsch PetscCall(PetscFree(auglag->Pscatter)); 332661095bbSAlp Dener if (tao->eq_constrained) { 3339566063dSJacob Faibussowitsch PetscCall(VecDestroy(&auglag->Y)); /* combo multipliers */ 3349566063dSJacob Faibussowitsch PetscCall(PetscFree(auglag->Yarr)); /* array of dual vectors */ 3359566063dSJacob Faibussowitsch PetscCall(VecDestroy(&auglag->C)); /* combo constraints */ 3369566063dSJacob Faibussowitsch PetscCall(ISDestroy(&auglag->Yis[0])); /* index set for Ye inside Y */ 3379566063dSJacob Faibussowitsch PetscCall(ISDestroy(&auglag->Yis[1])); /* index set for Yi inside Y */ 3389566063dSJacob Faibussowitsch PetscCall(PetscFree(auglag->Yis)); 3399566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&auglag->Yscatter[0])); 3409566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&auglag->Yscatter[1])); 3419566063dSJacob Faibussowitsch PetscCall(PetscFree(auglag->Yscatter)); 342661095bbSAlp Dener } 343661095bbSAlp Dener } 344661095bbSAlp Dener if (tao->bounded) { 3459566063dSJacob Faibussowitsch PetscCall(VecDestroy(&auglag->PL)); /* lower bounds for subsolver */ 3469566063dSJacob Faibussowitsch PetscCall(VecDestroy(&auglag->PU)); /* upper bounds for subsolver */ 347661095bbSAlp Dener } 348661095bbSAlp Dener } 3492e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoALMMGetType_C",NULL)); 3502e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoALMMSetType_C",NULL)); 3512e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoALMMGetSubsolver_C",NULL)); 3522e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoALMMSetSubsolver_C",NULL)); 3532e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoALMMGetMultipliers_C",NULL)); 3542e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoALMMSetMultipliers_C",NULL)); 3552e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoALMMGetPrimalIS_C",NULL)); 3562e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)tao,"TaoALMMGetDualIS_C",NULL)); 3579566063dSJacob Faibussowitsch PetscCall(PetscFree(tao->data)); 358661095bbSAlp Dener PetscFunctionReturn(0); 359661095bbSAlp Dener } 360661095bbSAlp Dener 361*dbbe0bcdSBarry Smith static PetscErrorCode TaoSetFromOptions_ALMM(Tao tao,PetscOptionItems *PetscOptionsObject) 362661095bbSAlp Dener { 363661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM*)tao->data; 364661095bbSAlp Dener PetscInt i; 365661095bbSAlp Dener 366661095bbSAlp Dener PetscFunctionBegin; 36771075aafSSatish Balay PetscOptionsHeadBegin(PetscOptionsObject,"Augmented Lagrangian multiplier method solves problems with general constraints by converting them into a sequence of unconstrained problems."); 3689566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_almm_mu_init","initial penalty parameter","",auglag->mu0,&auglag->mu0,NULL)); 3699566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_almm_mu_factor","increase factor for the penalty parameter","",auglag->mu_fac,&auglag->mu_fac,NULL)); 3709566063dSJacob 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)); 3719566063dSJacob 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)); 3729566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_almm_mu_max","maximum safeguard for penalty parameter updates","",auglag->mu_max,&auglag->mu_max,NULL)); 3739566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_almm_ye_min","minimum safeguard for equality multiplier updates","",auglag->ye_min,&auglag->ye_min,NULL)); 3749566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_almm_ye_max","maximum safeguard for equality multipliers updates","",auglag->ye_max,&auglag->ye_max,NULL)); 3759566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_almm_yi_min","minimum safeguard for inequality multipliers updates","",auglag->yi_min,&auglag->yi_min,NULL)); 3769566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tao_almm_yi_max","maximum safeguard for inequality multipliers updates","",auglag->yi_max,&auglag->yi_max,NULL)); 3779566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-tao_almm_type","augmented Lagrangian formulation type for the subproblem","TaoALMMType",TaoALMMTypes,(PetscEnum)auglag->type,(PetscEnum*)&auglag->type,NULL)); 378d0609cedSBarry Smith PetscOptionsHeadEnd(); 3799566063dSJacob Faibussowitsch PetscCall(TaoSetOptionsPrefix(auglag->subsolver,((PetscObject)tao)->prefix)); 3809566063dSJacob Faibussowitsch PetscCall(TaoAppendOptionsPrefix(auglag->subsolver,"tao_almm_subsolver_")); 3819566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(auglag->subsolver)); 382661095bbSAlp Dener for (i=0; i<tao->numbermonitors; i++) { 3839566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)tao->monitorcontext[i])); 3849566063dSJacob Faibussowitsch PetscCall(TaoSetMonitor(auglag->subsolver, tao->monitor[i], tao->monitorcontext[i], tao->monitordestroy[i])); 385661095bbSAlp Dener if (tao->monitor[i] == TaoMonitorDefault || tao->monitor[i] == TaoDefaultCMonitor || tao->monitor[i] == TaoDefaultGMonitor || tao->monitor[i] == TaoDefaultSMonitor) { 386661095bbSAlp Dener auglag->info = PETSC_TRUE; 387661095bbSAlp Dener } 388661095bbSAlp Dener } 389661095bbSAlp Dener PetscFunctionReturn(0); 390661095bbSAlp Dener } 391661095bbSAlp Dener 392661095bbSAlp Dener /* -------------------------------------------------------- */ 393661095bbSAlp Dener 394661095bbSAlp Dener /*MC 395661095bbSAlp Dener TaoALMM - Augmented Lagrangian multiplier method for solving nonlinear optimization problems with general constraints. 396661095bbSAlp Dener 397661095bbSAlp Dener Options Database Keys: 398661095bbSAlp Dener + -tao_almm_mu_init <real> - initial penalty parameter (default: 10.) 399661095bbSAlp Dener . -tao_almm_mu_factor <real> - increase factor for the penalty parameter (default: 100.) 400661095bbSAlp Dener . -tao_almm_mu_max <real> - maximum safeguard for penalty parameter updates (default: 1.e20) 401661095bbSAlp Dener . -tao_almm_mu_power_good <real> - exponential for penalty parameter when multiplier update is accepted (default: 0.9) 402661095bbSAlp Dener . -tao_almm_mu_power_bad <real> - exponential for penalty parameter when multiplier update is rejected (default: 0.1) 403661095bbSAlp Dener . -tao_almm_ye_min <real> - minimum safeguard for equality multiplier updates (default: -1.e20) 404661095bbSAlp Dener . -tao_almm_ye_max <real> - maximum safeguard for equality multiplier updates (default: 1.e20) 405661095bbSAlp Dener . -tao_almm_yi_min <real> - minimum safeguard for inequality multiplier updates (default: -1.e20) 406661095bbSAlp Dener . -tao_almm_yi_max <real> - maximum safeguard for inequality multiplier updates (default: 1.e20) 407661095bbSAlp Dener - -tao_almm_type <classic,phr> - change formulation of the augmented Lagrangian merit function for the subproblem (default: classic) 408661095bbSAlp Dener 409661095bbSAlp Dener Level: beginner 410661095bbSAlp Dener 411661095bbSAlp Dener Notes: 412661095bbSAlp Dener This method converts a constrained problem into a sequence of unconstrained problems via the augmented 413661095bbSAlp Dener Lagrangian merit function. Bound constraints are pushed down to the subproblem without any modifications. 414661095bbSAlp Dener 415661095bbSAlp Dener Two formulations are offered for the subproblem: canonical Hestenes-Powell augmented Lagrangian with slack 416661095bbSAlp Dener variables for inequality constraints, and a slack-less Powell-Hestenes-Rockafellar (PHR) formulation utilizing a 417661095bbSAlp Dener pointwise max() penalty on inequality constraints. The canonical augmented Lagrangian formulation typically 418661095bbSAlp Dener converges faster for most problems. However, PHR may be desirable for problems featuring a large number 419661095bbSAlp Dener of inequality constraints because it avoids inflating the size of the subproblem with slack variables. 420661095bbSAlp Dener 421661095bbSAlp Dener The subproblem is solved using a nested first-order TAO solver. The user can retrieve a pointer to 422661095bbSAlp Dener the subsolver via TaoALMMGetSubsolver() or pass command line arguments to it using the 423661095bbSAlp Dener "-tao_almm_subsolver_" prefix. Currently, TaoALMM does not support second-order methods for the 424661095bbSAlp Dener subproblem. It is also highly recommended that the subsolver chosen by the user utilize a trust-region 425661095bbSAlp Dener strategy for globalization (default: TAOBQNKTR) especially if the outer problem features bound constraints. 426661095bbSAlp Dener 427661095bbSAlp Dener .vb 428661095bbSAlp Dener while unconverged 429661095bbSAlp Dener solve argmin_x L(x) s.t. l <= x <= u 430661095bbSAlp Dener if ||c|| <= y_tol 431661095bbSAlp Dener if ||c|| <= c_tol && ||Lgrad|| <= g_tol: 432661095bbSAlp Dener problem converged, return solution 433661095bbSAlp Dener else 434661095bbSAlp Dener constraints sufficiently improved 435661095bbSAlp Dener update multipliers and tighten tolerances 436661095bbSAlp Dener endif 437661095bbSAlp Dener else 438661095bbSAlp Dener constraints did not improve 439661095bbSAlp Dener update penalty and loosen tolerances 440661095bbSAlp Dener endif 441661095bbSAlp Dener endwhile 442661095bbSAlp Dener .ve 443661095bbSAlp Dener 444db781477SPatrick Sanan .seealso: `TaoALMMGetType()`, `TaoALMMSetType()`, `TaoALMMSetSubsolver()`, `TaoALMMGetSubsolver()`, 445db781477SPatrick Sanan `TaoALMMGetMultipliers()`, `TaoALMMSetMultipliers()`, `TaoALMMGetPrimalIS()`, `TaoALMMGetDualIS()` 446661095bbSAlp Dener M*/ 447661095bbSAlp Dener PETSC_EXTERN PetscErrorCode TaoCreate_ALMM(Tao tao) 448661095bbSAlp Dener { 449661095bbSAlp Dener TAO_ALMM *auglag; 450661095bbSAlp Dener 451661095bbSAlp Dener PetscFunctionBegin; 4529566063dSJacob Faibussowitsch PetscCall(PetscNewLog(tao, &auglag)); 453661095bbSAlp Dener 454661095bbSAlp Dener tao->ops->destroy = TaoDestroy_ALMM; 455661095bbSAlp Dener tao->ops->setup = TaoSetUp_ALMM; 456661095bbSAlp Dener tao->ops->setfromoptions = TaoSetFromOptions_ALMM; 457661095bbSAlp Dener tao->ops->view = TaoView_ALMM; 458661095bbSAlp Dener tao->ops->solve = TaoSolve_ALMM; 459661095bbSAlp Dener 460661095bbSAlp Dener tao->gatol = 1.e-5; 461661095bbSAlp Dener tao->grtol = 0.0; 462661095bbSAlp Dener tao->gttol = 0.0; 463661095bbSAlp Dener tao->catol = 1.e-5; 464661095bbSAlp Dener tao->crtol = 0.0; 465661095bbSAlp Dener 466661095bbSAlp Dener tao->data = (void*)auglag; 467661095bbSAlp Dener auglag->parent = tao; 468661095bbSAlp Dener auglag->mu0 = 10.0; 469661095bbSAlp Dener auglag->mu = auglag->mu0; 470661095bbSAlp Dener auglag->mu_fac = 10.0; 471661095bbSAlp Dener auglag->mu_max = PETSC_INFINITY; 472661095bbSAlp Dener auglag->mu_pow_good = 0.9; 473661095bbSAlp Dener auglag->mu_pow_bad = 0.1; 474661095bbSAlp Dener auglag->ye_min = PETSC_NINFINITY; 475661095bbSAlp Dener auglag->ye_max = PETSC_INFINITY; 476661095bbSAlp Dener auglag->yi_min = PETSC_NINFINITY; 477661095bbSAlp Dener auglag->yi_max = PETSC_INFINITY; 478661095bbSAlp Dener auglag->ytol0 = 0.1/PetscPowReal(auglag->mu0, auglag->mu_pow_bad); 479661095bbSAlp Dener auglag->ytol = auglag->ytol0; 480661095bbSAlp Dener auglag->gtol0 = 1.0/auglag->mu0; 481661095bbSAlp Dener auglag->gtol = auglag->gtol0; 482661095bbSAlp Dener 483661095bbSAlp Dener auglag->sub_obj = TaoALMMComputeAugLagAndGradient_Private; 484661095bbSAlp Dener auglag->type = TAO_ALMM_CLASSIC; 485661095bbSAlp Dener auglag->info = PETSC_FALSE; 486661095bbSAlp Dener 4879566063dSJacob Faibussowitsch PetscCall(TaoCreate(PetscObjectComm((PetscObject)tao),&auglag->subsolver)); 4889566063dSJacob Faibussowitsch PetscCall(TaoSetType(auglag->subsolver, TAOBQNKTR)); 4899566063dSJacob Faibussowitsch PetscCall(TaoSetTolerances(auglag->subsolver, auglag->gtol, 0.0, 0.0)); 4909566063dSJacob Faibussowitsch PetscCall(TaoSetMaximumIterations(auglag->subsolver, 1000)); 4919566063dSJacob Faibussowitsch PetscCall(TaoSetMaximumFunctionEvaluations(auglag->subsolver, 10000)); 4929566063dSJacob Faibussowitsch PetscCall(TaoSetFunctionLowerBound(auglag->subsolver, PETSC_NINFINITY)); 4939566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)auglag->subsolver,(PetscObject)tao,1)); 494661095bbSAlp Dener 4959566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetType_C", TaoALMMGetType_Private)); 4969566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMSetType_C", TaoALMMSetType_Private)); 4979566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetSubsolver_C", TaoALMMGetSubsolver_Private)); 4989566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMSetSubsolver_C", TaoALMMSetSubsolver_Private)); 4999566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetMultipliers_C", TaoALMMGetMultipliers_Private)); 5009566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMSetMultipliers_C", TaoALMMSetMultipliers_Private)); 5019566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetPrimalIS_C", TaoALMMGetPrimalIS_Private)); 5029566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetDualIS_C", TaoALMMGetDualIS_Private)); 503661095bbSAlp Dener PetscFunctionReturn(0); 504661095bbSAlp Dener } 505661095bbSAlp Dener 506661095bbSAlp Dener static PetscErrorCode TaoALMMCombinePrimal_Private(Tao tao, Vec X, Vec S, Vec P) 507661095bbSAlp Dener { 508661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM*)tao->data; 509661095bbSAlp Dener 510661095bbSAlp Dener PetscFunctionBegin; 511661095bbSAlp Dener if (tao->ineq_constrained) { 5129566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(auglag->Pscatter[0], X, P, INSERT_VALUES, SCATTER_REVERSE)); 5139566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(auglag->Pscatter[0], X, P, INSERT_VALUES, SCATTER_REVERSE)); 5149566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(auglag->Pscatter[1], S, P, INSERT_VALUES, SCATTER_REVERSE)); 5159566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(auglag->Pscatter[1], S, P, INSERT_VALUES, SCATTER_REVERSE)); 516661095bbSAlp Dener } else { 5179566063dSJacob Faibussowitsch PetscCall(VecCopy(X, P)); 518661095bbSAlp Dener } 519661095bbSAlp Dener PetscFunctionReturn(0); 520661095bbSAlp Dener } 521661095bbSAlp Dener 522661095bbSAlp Dener static PetscErrorCode TaoALMMCombineDual_Private(Tao tao, Vec EQ, Vec IN, Vec Y) 523661095bbSAlp Dener { 524661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM*)tao->data; 525661095bbSAlp Dener 526661095bbSAlp Dener PetscFunctionBegin; 527661095bbSAlp Dener if (tao->eq_constrained) { 528661095bbSAlp Dener if (tao->ineq_constrained) { 5299566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(auglag->Yscatter[0], EQ, Y, INSERT_VALUES, SCATTER_REVERSE)); 5309566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(auglag->Yscatter[0], EQ, Y, INSERT_VALUES, SCATTER_REVERSE)); 5319566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(auglag->Yscatter[1], IN, Y, INSERT_VALUES, SCATTER_REVERSE)); 5329566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(auglag->Yscatter[1], IN, Y, INSERT_VALUES, SCATTER_REVERSE)); 533661095bbSAlp Dener } else { 5349566063dSJacob Faibussowitsch PetscCall(VecCopy(EQ, Y)); 535661095bbSAlp Dener } 536661095bbSAlp Dener } else { 5379566063dSJacob Faibussowitsch PetscCall(VecCopy(IN, Y)); 538661095bbSAlp Dener } 539661095bbSAlp Dener PetscFunctionReturn(0); 540661095bbSAlp Dener } 541661095bbSAlp Dener 542661095bbSAlp Dener static PetscErrorCode TaoALMMSplitPrimal_Private(Tao tao, Vec P, Vec X, Vec S) 543661095bbSAlp Dener { 544661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM*)tao->data; 545661095bbSAlp Dener 546661095bbSAlp Dener PetscFunctionBegin; 547661095bbSAlp Dener if (tao->ineq_constrained) { 5489566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(auglag->Pscatter[0], P, X, INSERT_VALUES, SCATTER_FORWARD)); 5499566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(auglag->Pscatter[0], P, X, INSERT_VALUES, SCATTER_FORWARD)); 5509566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(auglag->Pscatter[1], P, S, INSERT_VALUES, SCATTER_FORWARD)); 5519566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(auglag->Pscatter[1], P, S, INSERT_VALUES, SCATTER_FORWARD)); 552661095bbSAlp Dener } else { 5539566063dSJacob Faibussowitsch PetscCall(VecCopy(P, X)); 554661095bbSAlp Dener } 555661095bbSAlp Dener PetscFunctionReturn(0); 556661095bbSAlp Dener } 557661095bbSAlp Dener 558661095bbSAlp Dener /* this assumes that the latest constraints are stored in Ce and Ci, and also combined in C */ 559661095bbSAlp Dener static PetscErrorCode TaoALMMComputeOptimalityNorms_Private(Tao tao) 560661095bbSAlp Dener { 561661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM*)tao->data; 562661095bbSAlp Dener 563661095bbSAlp Dener PetscFunctionBegin; 564661095bbSAlp Dener /* if bounded, project the gradient */ 5651baa6e33SBarry Smith if (tao->bounded) PetscCall(VecBoundGradientProjection(auglag->LgradX, auglag->Px, tao->XL, tao->XU, auglag->LgradX)); 566661095bbSAlp Dener if (auglag->type == TAO_ALMM_PHR) { 5679566063dSJacob Faibussowitsch PetscCall(VecNorm(auglag->LgradX, NORM_INFINITY, &auglag->gnorm)); 568661095bbSAlp Dener auglag->cenorm = 0.0; 569661095bbSAlp Dener if (tao->eq_constrained) { 5709566063dSJacob Faibussowitsch PetscCall(VecNorm(auglag->Ce, NORM_INFINITY, &auglag->cenorm)); 571661095bbSAlp Dener } 572661095bbSAlp Dener auglag->cinorm = 0.0; 573661095bbSAlp Dener if (tao->ineq_constrained) { 5749566063dSJacob Faibussowitsch PetscCall(VecCopy(auglag->Yi, auglag->Ciwork)); 5759566063dSJacob Faibussowitsch PetscCall(VecScale(auglag->Ciwork, -1.0/auglag->mu)); 5769566063dSJacob Faibussowitsch PetscCall(VecPointwiseMax(auglag->Ciwork, auglag->Ci, auglag->Ciwork)); 5779566063dSJacob Faibussowitsch PetscCall(VecNorm(auglag->Ciwork, NORM_INFINITY, &auglag->cinorm)); 578661095bbSAlp Dener } 579661095bbSAlp Dener auglag->cnorm_old = auglag->cnorm; 580661095bbSAlp Dener auglag->cnorm = PetscMax(auglag->cenorm, auglag->cinorm); 581661095bbSAlp Dener auglag->ytol = auglag->ytol0 * auglag->cnorm_old; 582661095bbSAlp Dener } else { 5839566063dSJacob Faibussowitsch PetscCall(VecNorm(auglag->LgradX, NORM_2, &auglag->gnorm)); 5849566063dSJacob Faibussowitsch PetscCall(VecNorm(auglag->C, NORM_2, &auglag->cnorm)); 585661095bbSAlp Dener } 586661095bbSAlp Dener PetscFunctionReturn(0); 587661095bbSAlp Dener } 588661095bbSAlp Dener 589661095bbSAlp Dener static PetscErrorCode TaoALMMEvaluateIterate_Private(Tao tao, Vec P) 590661095bbSAlp Dener { 591661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM*)tao->data; 592661095bbSAlp Dener 593661095bbSAlp Dener PetscFunctionBegin; 594661095bbSAlp Dener /* split solution into primal and slack components */ 5959566063dSJacob Faibussowitsch PetscCall(TaoALMMSplitPrimal_Private(tao, auglag->P, auglag->Px, auglag->Ps)); 596661095bbSAlp Dener 597661095bbSAlp Dener /* compute f, df/dx and the constraints */ 5989566063dSJacob Faibussowitsch PetscCall(TaoComputeObjectiveAndGradient(tao, auglag->Px, &auglag->fval, auglag->LgradX)); 599661095bbSAlp Dener if (tao->eq_constrained) { 6009566063dSJacob Faibussowitsch PetscCall(TaoComputeEqualityConstraints(tao, auglag->Px, auglag->Ce)); 6019566063dSJacob Faibussowitsch PetscCall(TaoComputeJacobianEquality(tao, auglag->Px, auglag->Ae, auglag->Ae)); 602661095bbSAlp Dener } 603661095bbSAlp Dener if (tao->ineq_constrained) { 6049566063dSJacob Faibussowitsch PetscCall(TaoComputeInequalityConstraints(tao, auglag->Px, auglag->Ci)); 6059566063dSJacob Faibussowitsch PetscCall(TaoComputeJacobianInequality(tao, auglag->Px, auglag->Ai, auglag->Ai)); 606661095bbSAlp Dener switch (auglag->type) { 6077721a36fSStefano Zampini case TAO_ALMM_CLASSIC: 608661095bbSAlp Dener /* classic formulation converts inequality to equality constraints via slack variables */ 6099566063dSJacob Faibussowitsch PetscCall(VecAXPY(auglag->Ci, -1.0, auglag->Ps)); 610661095bbSAlp Dener break; 6117721a36fSStefano Zampini case TAO_ALMM_PHR: 612661095bbSAlp Dener /* PHR is based on Ci <= 0 while TAO defines Ci >= 0 so we hit it with a negative sign */ 6139566063dSJacob Faibussowitsch PetscCall(VecScale(auglag->Ci, -1.0)); 6149566063dSJacob Faibussowitsch PetscCall(MatScale(auglag->Ai, -1.0)); 615661095bbSAlp Dener break; 616661095bbSAlp Dener default: 617661095bbSAlp Dener break; 618661095bbSAlp Dener } 619661095bbSAlp Dener } 620661095bbSAlp Dener /* combine constraints into one vector */ 6219566063dSJacob Faibussowitsch PetscCall(TaoALMMCombineDual_Private(tao, auglag->Ce, auglag->Ci, auglag->C)); 622661095bbSAlp Dener PetscFunctionReturn(0); 623661095bbSAlp Dener } 624661095bbSAlp Dener 625661095bbSAlp Dener /* 626661095bbSAlp Dener Lphr = f + 0.5*mu*[ (Ce + Ye/mu)^T (Ce + Ye/mu) + pmin(0, Ci + Yi/mu)^T pmin(0, Ci + Yi/mu)] 627661095bbSAlp Dener 628661095bbSAlp Dener dLphr/dX = dF/dX + mu*[ (Ce + Ye/mu)^T Ae + pmin(0, Ci + Yi/mu)^T Ai] 629661095bbSAlp Dener 630661095bbSAlp Dener dLphr/dS = 0 631661095bbSAlp Dener */ 632661095bbSAlp Dener static PetscErrorCode TaoALMMComputePHRLagAndGradient_Private(Tao tao) 633661095bbSAlp Dener { 634661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM*)tao->data; 635661095bbSAlp Dener PetscReal eq_norm=0.0, ineq_norm=0.0; 636661095bbSAlp Dener 637661095bbSAlp Dener PetscFunctionBegin; 6389566063dSJacob Faibussowitsch PetscCall(TaoALMMEvaluateIterate_Private(tao, auglag->P)); 639661095bbSAlp Dener if (tao->eq_constrained) { 640661095bbSAlp Dener /* Ce_work = mu*(Ce + Ye/mu) */ 6419566063dSJacob Faibussowitsch PetscCall(VecWAXPY(auglag->Cework, 1.0/auglag->mu, auglag->Ye, auglag->Ce)); 6429566063dSJacob Faibussowitsch PetscCall(VecDot(auglag->Cework, auglag->Cework, &eq_norm)); /* contribution to scalar Lagrangian */ 6439566063dSJacob Faibussowitsch PetscCall(VecScale(auglag->Cework, auglag->mu)); 644661095bbSAlp Dener /* dL/dX += mu*(Ce + Ye/mu)^T Ae */ 6459566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(auglag->Ae, auglag->Cework, auglag->LgradX, auglag->LgradX)); 646661095bbSAlp Dener } 647661095bbSAlp Dener if (tao->ineq_constrained) { 648661095bbSAlp Dener /* Ci_work = mu * pmax(0, Ci + Yi/mu) where pmax() is pointwise max() */ 6499566063dSJacob Faibussowitsch PetscCall(VecWAXPY(auglag->Ciwork, 1.0/auglag->mu, auglag->Yi, auglag->Ci)); 6509566063dSJacob Faibussowitsch PetscCall(VecPointwiseMax(auglag->Ciwork, auglag->Cizero, auglag->Ciwork)); 6519566063dSJacob Faibussowitsch PetscCall(VecDot(auglag->Ciwork, auglag->Ciwork, &ineq_norm)); /* contribution to scalar Lagrangian */ 652661095bbSAlp Dener /* dL/dX += mu * pmax(0, Ci + Yi/mu)^T Ai */ 6539566063dSJacob Faibussowitsch PetscCall(VecScale(auglag->Ciwork, auglag->mu)); 6549566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(auglag->Ai, auglag->Ciwork, auglag->LgradX, auglag->LgradX)); 655a5b23f4aSJose E. Roman /* dL/dS = 0 because there are no slacks in PHR */ 6569566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(auglag->LgradS)); 657661095bbSAlp Dener } 658661095bbSAlp Dener /* combine gradient together */ 6599566063dSJacob Faibussowitsch PetscCall(TaoALMMCombinePrimal_Private(tao, auglag->LgradX, auglag->LgradS, auglag->G)); 660661095bbSAlp 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)] */ 6615f80ce2aSJacob Faibussowitsch auglag->Lval = auglag->fval + 0.5*auglag->mu*(eq_norm + ineq_norm); 662661095bbSAlp Dener PetscFunctionReturn(0); 663661095bbSAlp Dener } 664661095bbSAlp Dener 665661095bbSAlp Dener /* 666661095bbSAlp Dener Lc = F + Ye^TCe + Yi^T(Ci - S) + 0.5*mu*[Ce^TCe + (Ci - S)^T(Ci - S)] 667661095bbSAlp Dener 668661095bbSAlp Dener dLc/dX = dF/dX + Ye^TAe + Yi^TAi + 0.5*mu*[Ce^TAe + (Ci - S)^TAi] 669661095bbSAlp Dener 670661095bbSAlp Dener dLc/dS = -[Yi + mu*(Ci - S)] 671661095bbSAlp Dener */ 672661095bbSAlp Dener static PetscErrorCode TaoALMMComputeAugLagAndGradient_Private(Tao tao) 673661095bbSAlp Dener { 674661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM*)tao->data; 675661095bbSAlp Dener PetscReal yeTce=0.0, yiTcims=0.0, ceTce=0.0, cimsTcims=0.0; 676661095bbSAlp Dener 677661095bbSAlp Dener PetscFunctionBegin; 6789566063dSJacob Faibussowitsch PetscCall(TaoALMMEvaluateIterate_Private(tao, auglag->P)); 679661095bbSAlp Dener if (tao->eq_constrained) { 680661095bbSAlp Dener /* compute scalar contributions */ 6819566063dSJacob Faibussowitsch PetscCall(VecDot(auglag->Ye, auglag->Ce, &yeTce)); 6829566063dSJacob Faibussowitsch PetscCall(VecDot(auglag->Ce, auglag->Ce, &ceTce)); 683661095bbSAlp Dener /* dL/dX += ye^T Ae */ 6849566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(auglag->Ae, auglag->Ye, auglag->LgradX, auglag->LgradX)); 685661095bbSAlp Dener /* dL/dX += mu * ce^T Ae */ 6869566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(auglag->Ae, auglag->Ce, auglag->Xwork)); 6879566063dSJacob Faibussowitsch PetscCall(VecAXPY(auglag->LgradX, auglag->mu, auglag->Xwork)); 688661095bbSAlp Dener } 689661095bbSAlp Dener if (tao->ineq_constrained) { 690661095bbSAlp Dener /* compute scalar contributions */ 6919566063dSJacob Faibussowitsch PetscCall(VecDot(auglag->Yi, auglag->Ci, &yiTcims)); 6929566063dSJacob Faibussowitsch PetscCall(VecDot(auglag->Ci, auglag->Ci, &cimsTcims)); 693661095bbSAlp Dener /* dL/dX += yi^T Ai */ 6949566063dSJacob Faibussowitsch PetscCall(MatMultTransposeAdd(auglag->Ai, auglag->Yi, auglag->LgradX, auglag->LgradX)); 695661095bbSAlp Dener /* dL/dX += mu * (ci - s)^T Ai */ 6969566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(auglag->Ai, auglag->Ci, auglag->Xwork)); 6979566063dSJacob Faibussowitsch PetscCall(VecAXPY(auglag->LgradX, auglag->mu, auglag->Xwork)); 698661095bbSAlp Dener /* dL/dS = -[yi + mu*(ci - s)] */ 6999566063dSJacob Faibussowitsch PetscCall(VecWAXPY(auglag->LgradS, auglag->mu, auglag->Ci, auglag->Yi)); 7009566063dSJacob Faibussowitsch PetscCall(VecScale(auglag->LgradS, -1.0)); 701661095bbSAlp Dener } 702661095bbSAlp Dener /* combine gradient together */ 7039566063dSJacob Faibussowitsch PetscCall(TaoALMMCombinePrimal_Private(tao, auglag->LgradX, auglag->LgradS, auglag->G)); 704661095bbSAlp Dener /* compute L = f + ye^T ce + yi^T (ci - s) + 0.5*mu*||ce||^2 + 0.5*mu*||ci - s||^2 */ 705661095bbSAlp Dener auglag->Lval = auglag->fval + yeTce + yiTcims + 0.5*auglag->mu*(ceTce + cimsTcims); 706661095bbSAlp Dener PetscFunctionReturn(0); 707661095bbSAlp Dener } 708661095bbSAlp Dener 7097721a36fSStefano Zampini PetscErrorCode TaoALMMSubsolverObjective_Private(Tao tao, Vec P, PetscReal *Lval, void *ctx) 7107721a36fSStefano Zampini { 7117721a36fSStefano Zampini TAO_ALMM *auglag = (TAO_ALMM*)ctx; 7127721a36fSStefano Zampini 7137721a36fSStefano Zampini PetscFunctionBegin; 7149566063dSJacob Faibussowitsch PetscCall(VecCopy(P, auglag->P)); 7159566063dSJacob Faibussowitsch PetscCall((*auglag->sub_obj)(auglag->parent)); 7167721a36fSStefano Zampini *Lval = auglag->Lval; 7177721a36fSStefano Zampini PetscFunctionReturn(0); 7187721a36fSStefano Zampini } 7197721a36fSStefano Zampini 720661095bbSAlp Dener PetscErrorCode TaoALMMSubsolverObjectiveAndGradient_Private(Tao tao, Vec P, PetscReal *Lval, Vec G, void *ctx) 721661095bbSAlp Dener { 722661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM*)ctx; 723661095bbSAlp Dener 724661095bbSAlp Dener PetscFunctionBegin; 7259566063dSJacob Faibussowitsch PetscCall(VecCopy(P, auglag->P)); 7269566063dSJacob Faibussowitsch PetscCall((*auglag->sub_obj)(auglag->parent)); 7279566063dSJacob Faibussowitsch PetscCall(VecCopy(auglag->G, G)); 728661095bbSAlp Dener *Lval = auglag->Lval; 729661095bbSAlp Dener PetscFunctionReturn(0); 730661095bbSAlp Dener } 731