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 PetscErrorCode ierr; 19661095bbSAlp Dener 20661095bbSAlp Dener PetscFunctionBegin; 21661095bbSAlp Dener /* reset initial multiplier/slack guess */ 22661095bbSAlp Dener if (!tao->recycle) { 23661095bbSAlp Dener if (tao->ineq_constrained) { 24661095bbSAlp Dener ierr = VecZeroEntries(auglag->Ps);CHKERRQ(ierr); 25661095bbSAlp Dener ierr = TaoALMMCombinePrimal_Private(tao, auglag->Px, auglag->Ps, auglag->P);CHKERRQ(ierr); 26661095bbSAlp Dener ierr = VecZeroEntries(auglag->Yi);CHKERRQ(ierr); 27661095bbSAlp Dener } 28661095bbSAlp Dener if (tao->eq_constrained) { 29661095bbSAlp Dener ierr = VecZeroEntries(auglag->Ye);CHKERRQ(ierr); 30661095bbSAlp Dener } 31661095bbSAlp Dener } 32661095bbSAlp Dener 33661095bbSAlp Dener /* compute initial nonlinear Lagrangian and its derivatives */ 34661095bbSAlp Dener ierr = (*auglag->sub_obj)(tao);CHKERRQ(ierr); 35661095bbSAlp Dener ierr = TaoALMMComputeOptimalityNorms_Private(tao);CHKERRQ(ierr); 36661095bbSAlp Dener /* print initial step and check convergence */ 377d3de750SJacob Faibussowitsch ierr = PetscInfo(tao,"Solving with %s formulation\n",TaoALMMTypes[auglag->type]);CHKERRQ(ierr); 38661095bbSAlp Dener ierr = TaoLogConvergenceHistory(tao, auglag->Lval, auglag->gnorm, auglag->cnorm, tao->ksp_its);CHKERRQ(ierr); 39661095bbSAlp Dener ierr = TaoMonitor(tao, tao->niter, auglag->fval, auglag->gnorm, auglag->cnorm, 0.0);CHKERRQ(ierr); 40661095bbSAlp Dener ierr = (*tao->ops->convergencetest)(tao, tao->cnvP);CHKERRQ(ierr); 41661095bbSAlp Dener /* set initial penalty factor and inner solver tolerance */ 42661095bbSAlp Dener switch (auglag->type) { 437721a36fSStefano Zampini case TAO_ALMM_CLASSIC: 44661095bbSAlp Dener auglag->mu = auglag->mu0; 45661095bbSAlp Dener break; 467721a36fSStefano Zampini case TAO_ALMM_PHR: 47661095bbSAlp Dener auglag->cenorm = 0.0; 48661095bbSAlp Dener if (tao->eq_constrained) { 49661095bbSAlp Dener ierr = VecDot(auglag->Ce, auglag->Ce, &auglag->cenorm);CHKERRQ(ierr); 50661095bbSAlp Dener } 51661095bbSAlp Dener auglag->cinorm = 0.0; 52661095bbSAlp Dener if (tao->ineq_constrained) { 53661095bbSAlp Dener ierr = VecCopy(auglag->Ci, auglag->Ciwork);CHKERRQ(ierr); 54661095bbSAlp Dener ierr = VecScale(auglag->Ciwork, -1.0);CHKERRQ(ierr); 55661095bbSAlp Dener ierr = VecPointwiseMax(auglag->Ciwork, auglag->Cizero, auglag->Ciwork);CHKERRQ(ierr); 56661095bbSAlp Dener ierr = VecDot(auglag->Ciwork, auglag->Ciwork, &auglag->cinorm);CHKERRQ(ierr); 57661095bbSAlp Dener } 58661095bbSAlp Dener /* determine initial penalty factor based on the balance of constraint violation and objective function value */ 59661095bbSAlp Dener auglag->mu = PetscMax(1.e-6, PetscMin(10.0, 2.0*PetscAbsReal(auglag->fval)/(auglag->cenorm + auglag->cinorm))); 60661095bbSAlp Dener break; 61661095bbSAlp Dener default: 62661095bbSAlp Dener break; 63661095bbSAlp Dener } 64661095bbSAlp Dener auglag->gtol = auglag->gtol0; 657d3de750SJacob Faibussowitsch ierr = PetscInfo(tao,"Initial penalty: %.2f\n",auglag->mu);CHKERRQ(ierr); 66661095bbSAlp Dener 67661095bbSAlp Dener /* start aug-lag outer loop */ 68661095bbSAlp Dener while (tao->reason == TAO_CONTINUE_ITERATING) { 69661095bbSAlp Dener ++tao->niter; 70661095bbSAlp Dener /* update subsolver tolerance */ 717d3de750SJacob Faibussowitsch ierr = PetscInfo(tao,"Subsolver tolerance: ||G|| <= %e\n",auglag->gtol);CHKERRQ(ierr); 72661095bbSAlp Dener ierr = TaoSetTolerances(auglag->subsolver, auglag->gtol, 0.0, 0.0);CHKERRQ(ierr); 73661095bbSAlp Dener /* solve the bound-constrained or unconstrained subproblem */ 74661095bbSAlp Dener ierr = TaoSolve(auglag->subsolver);CHKERRQ(ierr); 75661095bbSAlp Dener ierr = TaoGetConvergedReason(auglag->subsolver, &reason);CHKERRQ(ierr); 76661095bbSAlp Dener tao->ksp_its += auglag->subsolver->ksp_its; 77661095bbSAlp Dener if (reason != TAO_CONVERGED_GATOL) { 787d3de750SJacob Faibussowitsch ierr = PetscInfo(tao,"Subsolver failed to converge, reason: %s\n",TaoConvergedReasons[reason]);CHKERRQ(ierr); 79661095bbSAlp Dener } 80661095bbSAlp Dener /* evaluate solution and test convergence */ 81661095bbSAlp Dener ierr = (*auglag->sub_obj)(tao);CHKERRQ(ierr); 82661095bbSAlp Dener ierr = TaoALMMComputeOptimalityNorms_Private(tao);CHKERRQ(ierr); 83661095bbSAlp Dener /* decide whether to update multipliers or not */ 84661095bbSAlp Dener updated = 0.0; 85661095bbSAlp Dener if (auglag->cnorm <= auglag->ytol) { 867d3de750SJacob Faibussowitsch ierr = PetscInfo(tao,"Multipliers updated: ||C|| <= %e\n",auglag->ytol);CHKERRQ(ierr); 87661095bbSAlp Dener /* constraints are good, update multipliers and convergence tolerances */ 88661095bbSAlp Dener if (tao->eq_constrained) { 89661095bbSAlp Dener ierr = VecAXPY(auglag->Ye, auglag->mu, auglag->Ce);CHKERRQ(ierr); 90661095bbSAlp Dener ierr = VecSet(auglag->Cework, auglag->ye_max);CHKERRQ(ierr); 91661095bbSAlp Dener ierr = VecPointwiseMin(auglag->Ye, auglag->Cework, auglag->Ye);CHKERRQ(ierr); 92661095bbSAlp Dener ierr = VecSet(auglag->Cework, auglag->ye_min);CHKERRQ(ierr); 93661095bbSAlp Dener ierr = VecPointwiseMax(auglag->Ye, auglag->Cework, auglag->Ye);CHKERRQ(ierr); 94661095bbSAlp Dener } 95661095bbSAlp Dener if (tao->ineq_constrained) { 96661095bbSAlp Dener ierr = VecAXPY(auglag->Yi, auglag->mu, auglag->Ci);CHKERRQ(ierr); 97661095bbSAlp Dener ierr = VecSet(auglag->Ciwork, auglag->yi_max);CHKERRQ(ierr); 98661095bbSAlp Dener ierr = VecPointwiseMin(auglag->Yi, auglag->Ciwork, auglag->Yi);CHKERRQ(ierr); 99661095bbSAlp Dener ierr = VecSet(auglag->Ciwork, auglag->yi_min);CHKERRQ(ierr); 100661095bbSAlp Dener ierr = VecPointwiseMax(auglag->Yi, auglag->Ciwork, auglag->Yi);CHKERRQ(ierr); 101661095bbSAlp Dener } 102661095bbSAlp Dener /* tolerances are updated only for non-PHR methods */ 103661095bbSAlp Dener if (auglag->type != TAO_ALMM_PHR) { 104661095bbSAlp Dener auglag->ytol = PetscMax(tao->catol, auglag->ytol/PetscPowReal(auglag->mu, auglag->mu_pow_good)); 105661095bbSAlp Dener auglag->gtol = PetscMax(tao->gatol, auglag->gtol/auglag->mu); 106661095bbSAlp Dener } 107661095bbSAlp Dener updated = 1.0; 108661095bbSAlp Dener } else { 109661095bbSAlp Dener /* constraints are bad, update penalty factor */ 110661095bbSAlp Dener auglag->mu = PetscMin(auglag->mu_max, auglag->mu_fac*auglag->mu); 111661095bbSAlp Dener /* tolerances are reset only for non-PHR methods */ 112661095bbSAlp Dener if (auglag->type != TAO_ALMM_PHR) { 113661095bbSAlp Dener auglag->ytol = PetscMax(tao->catol, 0.1/PetscPowReal(auglag->mu, auglag->mu_pow_bad)); 114661095bbSAlp Dener auglag->gtol = PetscMax(tao->gatol, 1.0/auglag->mu); 115661095bbSAlp Dener } 1167d3de750SJacob Faibussowitsch ierr = PetscInfo(tao,"Penalty increased: mu = %.2f\n",auglag->mu);CHKERRQ(ierr); 117661095bbSAlp Dener } 118661095bbSAlp Dener ierr = TaoLogConvergenceHistory(tao, auglag->fval, auglag->gnorm, auglag->cnorm, tao->ksp_its);CHKERRQ(ierr); 119661095bbSAlp Dener ierr = TaoMonitor(tao, tao->niter, auglag->fval, auglag->gnorm, auglag->cnorm, updated);CHKERRQ(ierr); 120661095bbSAlp Dener ierr = (*tao->ops->convergencetest)(tao, tao->cnvP);CHKERRQ(ierr); 121661095bbSAlp Dener } 122661095bbSAlp Dener 123661095bbSAlp Dener PetscFunctionReturn(0); 124661095bbSAlp Dener } 125661095bbSAlp Dener 126661095bbSAlp Dener static PetscErrorCode TaoView_ALMM(Tao tao,PetscViewer viewer) 127661095bbSAlp Dener { 128661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM*)tao->data; 129661095bbSAlp Dener PetscBool isascii; 130661095bbSAlp Dener PetscErrorCode ierr; 131661095bbSAlp Dener 132661095bbSAlp Dener PetscFunctionBegin; 133661095bbSAlp Dener ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 134661095bbSAlp Dener ierr = TaoView(auglag->subsolver,viewer);CHKERRQ(ierr); 135661095bbSAlp Dener ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 136661095bbSAlp Dener ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr); 137661095bbSAlp Dener if (isascii) { 138661095bbSAlp Dener ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 139661095bbSAlp Dener ierr = PetscViewerASCIIPrintf(viewer, "ALMM Formulation Type: %s\n", TaoALMMTypes[auglag->type]);CHKERRQ(ierr); 140661095bbSAlp Dener ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 141661095bbSAlp Dener } 142661095bbSAlp Dener PetscFunctionReturn(0); 143661095bbSAlp Dener } 144661095bbSAlp Dener 145661095bbSAlp Dener static PetscErrorCode TaoSetUp_ALMM(Tao tao) 146661095bbSAlp Dener { 147661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM*)tao->data; 148661095bbSAlp Dener VecType vec_type; 149661095bbSAlp Dener Vec SL, SU; 150661095bbSAlp Dener PetscBool is_cg = PETSC_FALSE, is_lmvm = PETSC_FALSE; 151661095bbSAlp Dener PetscErrorCode ierr; 152661095bbSAlp Dener 153661095bbSAlp Dener PetscFunctionBegin; 154*3c859ba3SBarry 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."); 155*3c859ba3SBarry 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."); 156661095bbSAlp Dener ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr); 157661095bbSAlp Dener /* alias base vectors and create extras */ 158661095bbSAlp Dener ierr = VecGetType(tao->solution, &vec_type);CHKERRQ(ierr); 159661095bbSAlp Dener auglag->Px = tao->solution; 160661095bbSAlp Dener if (!tao->gradient) { /* base gradient */ 161661095bbSAlp Dener ierr = VecDuplicate(tao->solution, &tao->gradient);CHKERRQ(ierr); 162661095bbSAlp Dener } 163661095bbSAlp Dener auglag->LgradX = tao->gradient; 164661095bbSAlp Dener if (!auglag->Xwork) { /* opt var work vector */ 165661095bbSAlp Dener ierr = VecDuplicate(tao->solution, &auglag->Xwork);CHKERRQ(ierr); 166661095bbSAlp Dener } 167661095bbSAlp Dener if (tao->eq_constrained) { 168661095bbSAlp Dener auglag->Ce = tao->constraints_equality; 169661095bbSAlp Dener auglag->Ae = tao->jacobian_equality; 170661095bbSAlp Dener if (!auglag->Ye) { /* equality multipliers */ 171661095bbSAlp Dener ierr = VecDuplicate(auglag->Ce, &auglag->Ye);CHKERRQ(ierr); 172661095bbSAlp Dener } 173661095bbSAlp Dener if (!auglag->Cework) { 174661095bbSAlp Dener ierr = VecDuplicate(auglag->Ce, &auglag->Cework);CHKERRQ(ierr); 175661095bbSAlp Dener } 176661095bbSAlp Dener } 177661095bbSAlp Dener if (tao->ineq_constrained) { 178661095bbSAlp Dener auglag->Ci = tao->constraints_inequality; 179661095bbSAlp Dener auglag->Ai = tao->jacobian_inequality; 180661095bbSAlp Dener if (!auglag->Yi) { /* inequality multipliers */ 181661095bbSAlp Dener ierr = VecDuplicate(auglag->Ci, &auglag->Yi);CHKERRQ(ierr); 182661095bbSAlp Dener } 183661095bbSAlp Dener if (!auglag->Ciwork) { 184661095bbSAlp Dener ierr = VecDuplicate(auglag->Ci, &auglag->Ciwork);CHKERRQ(ierr); 185661095bbSAlp Dener } 186661095bbSAlp Dener if (!auglag->Cizero) { 187661095bbSAlp Dener ierr = VecDuplicate(auglag->Ci, &auglag->Cizero);CHKERRQ(ierr); 188661095bbSAlp Dener ierr = VecZeroEntries(auglag->Cizero);CHKERRQ(ierr); 189661095bbSAlp Dener } 190661095bbSAlp Dener if (!auglag->Ps) { /* slack vars */ 191661095bbSAlp Dener ierr = VecDuplicate(auglag->Ci, &auglag->Ps);CHKERRQ(ierr); 192661095bbSAlp Dener } 193661095bbSAlp Dener if (!auglag->LgradS) { /* slack component of Lagrangian gradient */ 194661095bbSAlp Dener ierr = VecDuplicate(auglag->Ci, &auglag->LgradS);CHKERRQ(ierr); 195661095bbSAlp Dener } 196661095bbSAlp Dener /* create vector for combined primal space and the associated communication objects */ 197661095bbSAlp Dener if (!auglag->P) { 198661095bbSAlp Dener ierr = PetscMalloc1(2, &auglag->Parr);CHKERRQ(ierr); 199661095bbSAlp Dener auglag->Parr[0] = auglag->Px; auglag->Parr[1] = auglag->Ps; 2001e1ea65dSPierre Jolivet ierr = VecConcatenate(2, auglag->Parr, &auglag->P, &auglag->Pis);CHKERRQ(ierr); 201661095bbSAlp Dener ierr = PetscMalloc1(2, &auglag->Pscatter);CHKERRQ(ierr); 2021e1ea65dSPierre Jolivet ierr = VecScatterCreate(auglag->P, auglag->Pis[0], auglag->Px, NULL, &auglag->Pscatter[0]);CHKERRQ(ierr); 2031e1ea65dSPierre Jolivet ierr = VecScatterCreate(auglag->P, auglag->Pis[1], auglag->Ps, NULL, &auglag->Pscatter[1]);CHKERRQ(ierr); 204661095bbSAlp Dener } 205661095bbSAlp Dener if (tao->eq_constrained) { 206661095bbSAlp Dener /* create vector for combined dual space and the associated communication objects */ 207661095bbSAlp Dener if (!auglag->Y) { 208661095bbSAlp Dener ierr = PetscMalloc1(2, &auglag->Yarr);CHKERRQ(ierr); 209661095bbSAlp Dener auglag->Yarr[0] = auglag->Ye; auglag->Yarr[1] = auglag->Yi; 2101e1ea65dSPierre Jolivet ierr = VecConcatenate(2, auglag->Yarr, &auglag->Y, &auglag->Yis);CHKERRQ(ierr); 211661095bbSAlp Dener ierr = PetscMalloc1(2, &auglag->Yscatter);CHKERRQ(ierr); 2121e1ea65dSPierre Jolivet ierr = VecScatterCreate(auglag->Y, auglag->Yis[0], auglag->Ye, NULL, &auglag->Yscatter[0]);CHKERRQ(ierr); 2131e1ea65dSPierre Jolivet ierr = VecScatterCreate(auglag->Y, auglag->Yis[1], auglag->Yi, NULL, &auglag->Yscatter[1]);CHKERRQ(ierr); 214661095bbSAlp Dener } 215661095bbSAlp Dener if (!auglag->C) { 216661095bbSAlp Dener ierr = VecDuplicate(auglag->Y, &auglag->C);CHKERRQ(ierr); 217661095bbSAlp Dener } 218661095bbSAlp Dener } else { 219661095bbSAlp Dener if (!auglag->C) { 220661095bbSAlp Dener auglag->C = auglag->Ci; 221661095bbSAlp Dener } 222661095bbSAlp Dener if (!auglag->Y) { 223661095bbSAlp Dener auglag->Y = auglag->Yi; 224661095bbSAlp Dener } 225661095bbSAlp Dener } 226661095bbSAlp Dener } else { 227661095bbSAlp Dener if (!auglag->P) { 228661095bbSAlp Dener auglag->P = auglag->Px; 229661095bbSAlp Dener } 230661095bbSAlp Dener if (!auglag->G) { 231661095bbSAlp Dener auglag->G = auglag->LgradX; 232661095bbSAlp Dener } 233661095bbSAlp Dener if (!auglag->C) { 234661095bbSAlp Dener auglag->C = auglag->Ce; 235661095bbSAlp Dener } 236661095bbSAlp Dener if (!auglag->Y) { 237661095bbSAlp Dener auglag->Y = auglag->Ye; 238661095bbSAlp Dener } 239661095bbSAlp Dener } 240661095bbSAlp Dener /* initialize parameters */ 241661095bbSAlp Dener if (auglag->type == TAO_ALMM_PHR) { 242661095bbSAlp Dener auglag->mu_fac = 10.0; 243661095bbSAlp Dener auglag->yi_min = 0.0; 244661095bbSAlp Dener auglag->ytol0 = 0.5; 245661095bbSAlp Dener auglag->gtol0 = tao->gatol; 246661095bbSAlp Dener if (tao->gatol_changed && tao->catol_changed) { 247661095bbSAlp Dener ierr = PetscInfo(tao,"TAOALMM with PHR: different gradient and constraint tolerances are not supported, setting catol = gatol\n");CHKERRQ(ierr); 248661095bbSAlp Dener tao->catol = tao->gatol; 249661095bbSAlp Dener } 250661095bbSAlp Dener } 251661095bbSAlp Dener /* set the Lagrangian formulation type for the subsolver */ 252661095bbSAlp Dener switch (auglag->type) { 2537721a36fSStefano Zampini case TAO_ALMM_CLASSIC: 254661095bbSAlp Dener auglag->sub_obj = TaoALMMComputeAugLagAndGradient_Private; 255661095bbSAlp Dener break; 2567721a36fSStefano Zampini case TAO_ALMM_PHR: 257661095bbSAlp Dener auglag->sub_obj = TaoALMMComputePHRLagAndGradient_Private; 258661095bbSAlp Dener break; 259661095bbSAlp Dener default: 260661095bbSAlp Dener break; 261661095bbSAlp Dener } 262661095bbSAlp Dener /* set up the subsolver */ 263a82e8c82SStefano Zampini ierr = TaoSetSolution(auglag->subsolver, auglag->P);CHKERRQ(ierr); 264a82e8c82SStefano Zampini ierr = TaoSetObjective(auglag->subsolver, TaoALMMSubsolverObjective_Private, (void*)auglag);CHKERRQ(ierr); 265a82e8c82SStefano Zampini ierr = TaoSetObjectiveAndGradient(auglag->subsolver, NULL, TaoALMMSubsolverObjectiveAndGradient_Private, (void*)auglag);CHKERRQ(ierr); 266661095bbSAlp Dener if (tao->bounded) { 267661095bbSAlp Dener /* make sure that the subsolver is a bound-constrained method */ 268661095bbSAlp Dener ierr = PetscObjectTypeCompare((PetscObject)auglag->subsolver, TAOCG, &is_cg);CHKERRQ(ierr); 269661095bbSAlp Dener ierr = PetscObjectTypeCompare((PetscObject)auglag->subsolver, TAOLMVM, &is_lmvm);CHKERRQ(ierr); 270661095bbSAlp Dener if (is_cg) { 271661095bbSAlp Dener ierr = TaoSetType(auglag->subsolver, TAOBNCG);CHKERRQ(ierr); 272661095bbSAlp Dener ierr = PetscInfo(tao,"TAOCG detected for bound-constrained problem, switching to TAOBNCG instead.");CHKERRQ(ierr); 273661095bbSAlp Dener } 274661095bbSAlp Dener if (is_lmvm) { 275661095bbSAlp Dener ierr = TaoSetType(auglag->subsolver, TAOBQNLS);CHKERRQ(ierr); 276661095bbSAlp Dener ierr = PetscInfo(tao,"TAOLMVM detected for bound-constrained problem, switching to TAOBQNLS instead.");CHKERRQ(ierr); 277661095bbSAlp Dener } 278661095bbSAlp Dener /* create lower and upper bound clone vectors for subsolver */ 279661095bbSAlp Dener if (!auglag->PL) { 280661095bbSAlp Dener ierr = VecDuplicate(auglag->P, &auglag->PL);CHKERRQ(ierr); 281661095bbSAlp Dener } 282661095bbSAlp Dener if (!auglag->PU) { 283661095bbSAlp Dener ierr = VecDuplicate(auglag->P, &auglag->PU);CHKERRQ(ierr); 284661095bbSAlp Dener } 285661095bbSAlp Dener if (tao->ineq_constrained) { 286661095bbSAlp Dener /* create lower and upper bounds for slack, set lower to 0 */ 287661095bbSAlp Dener ierr = VecDuplicate(auglag->Ci, &SL);CHKERRQ(ierr); 288661095bbSAlp Dener ierr = VecSet(SL, 0.0);CHKERRQ(ierr); 289661095bbSAlp Dener ierr = VecDuplicate(auglag->Ci, &SU);CHKERRQ(ierr); 290661095bbSAlp Dener ierr = VecSet(SU, PETSC_INFINITY);CHKERRQ(ierr); 291661095bbSAlp Dener /* combine opt var bounds with slack bounds */ 292661095bbSAlp Dener ierr = TaoALMMCombinePrimal_Private(tao, tao->XL, SL, auglag->PL);CHKERRQ(ierr); 293661095bbSAlp Dener ierr = TaoALMMCombinePrimal_Private(tao, tao->XU, SU, auglag->PU);CHKERRQ(ierr); 294661095bbSAlp Dener /* destroy work vectors */ 295661095bbSAlp Dener ierr = VecDestroy(&SL);CHKERRQ(ierr); 296661095bbSAlp Dener ierr = VecDestroy(&SU);CHKERRQ(ierr); 297661095bbSAlp Dener } else { 298661095bbSAlp Dener /* no inequality constraints, just copy bounds into the subsolver */ 299661095bbSAlp Dener ierr = VecCopy(tao->XL, auglag->PL);CHKERRQ(ierr); 300661095bbSAlp Dener ierr = VecCopy(tao->XU, auglag->PU);CHKERRQ(ierr); 301661095bbSAlp Dener } 302661095bbSAlp Dener ierr = TaoSetVariableBounds(auglag->subsolver, auglag->PL, auglag->PU);CHKERRQ(ierr); 303661095bbSAlp Dener } 304661095bbSAlp Dener ierr = TaoSetUp(auglag->subsolver);CHKERRQ(ierr); 305661095bbSAlp Dener auglag->G = auglag->subsolver->gradient; 306661095bbSAlp Dener 307661095bbSAlp Dener PetscFunctionReturn(0); 308661095bbSAlp Dener } 309661095bbSAlp Dener 310661095bbSAlp Dener static PetscErrorCode TaoDestroy_ALMM(Tao tao) 311661095bbSAlp Dener { 312661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM*)tao->data; 313661095bbSAlp Dener PetscErrorCode ierr; 314661095bbSAlp Dener 315661095bbSAlp Dener PetscFunctionBegin; 316661095bbSAlp Dener ierr = TaoDestroy(&auglag->subsolver);CHKERRQ(ierr); 317661095bbSAlp Dener if (tao->setupcalled) { 318661095bbSAlp Dener ierr = VecDestroy(&auglag->Xwork);CHKERRQ(ierr); /* opt work */ 319661095bbSAlp Dener if (tao->eq_constrained) { 320661095bbSAlp Dener ierr = VecDestroy(&auglag->Ye);CHKERRQ(ierr); /* equality multipliers */ 321661095bbSAlp Dener ierr = VecDestroy(&auglag->Cework);CHKERRQ(ierr); /* equality work vector */ 322661095bbSAlp Dener } 323661095bbSAlp Dener if (tao->ineq_constrained) { 324661095bbSAlp Dener ierr = VecDestroy(&auglag->Ps);CHKERRQ(ierr); /* slack vars */ 325661095bbSAlp Dener auglag->Parr[0] = NULL; /* clear pointer to tao->solution, will be destroyed by TaoDestroy() shell */ 326661095bbSAlp Dener ierr = PetscFree(auglag->Parr);CHKERRQ(ierr); /* array of primal vectors */ 327661095bbSAlp Dener ierr = VecDestroy(&auglag->LgradS);CHKERRQ(ierr); /* slack grad */ 328661095bbSAlp Dener ierr = VecDestroy(&auglag->Cizero);CHKERRQ(ierr); /* zero vector for pointwise max */ 329661095bbSAlp Dener ierr = VecDestroy(&auglag->Yi);CHKERRQ(ierr); /* inequality multipliers */ 330661095bbSAlp Dener ierr = VecDestroy(&auglag->Ciwork);CHKERRQ(ierr); /* inequality work vector */ 331661095bbSAlp Dener ierr = VecDestroy(&auglag->P);CHKERRQ(ierr); /* combo primal */ 332661095bbSAlp Dener ierr = ISDestroy(&auglag->Pis[0]);CHKERRQ(ierr); /* index set for X inside P */ 333661095bbSAlp Dener ierr = ISDestroy(&auglag->Pis[1]);CHKERRQ(ierr); /* index set for S inside P */ 334661095bbSAlp Dener ierr = PetscFree(auglag->Pis);CHKERRQ(ierr); /* array of P index sets */ 335661095bbSAlp Dener ierr = VecScatterDestroy(&auglag->Pscatter[0]);CHKERRQ(ierr); 336661095bbSAlp Dener ierr = VecScatterDestroy(&auglag->Pscatter[1]);CHKERRQ(ierr); 337661095bbSAlp Dener ierr = PetscFree(auglag->Pscatter);CHKERRQ(ierr); 338661095bbSAlp Dener if (tao->eq_constrained) { 339661095bbSAlp Dener ierr = VecDestroy(&auglag->Y);CHKERRQ(ierr); /* combo multipliers */ 340661095bbSAlp Dener ierr = PetscFree(auglag->Yarr);CHKERRQ(ierr); /* array of dual vectors */ 341661095bbSAlp Dener ierr = VecDestroy(&auglag->C);CHKERRQ(ierr); /* combo constraints */ 342661095bbSAlp Dener ierr = ISDestroy(&auglag->Yis[0]);CHKERRQ(ierr); /* index set for Ye inside Y */ 343661095bbSAlp Dener ierr = ISDestroy(&auglag->Yis[1]);CHKERRQ(ierr); /* index set for Yi inside Y */ 344661095bbSAlp Dener ierr = PetscFree(auglag->Yis);CHKERRQ(ierr); 345661095bbSAlp Dener ierr = VecScatterDestroy(&auglag->Yscatter[0]);CHKERRQ(ierr); 346661095bbSAlp Dener ierr = VecScatterDestroy(&auglag->Yscatter[1]);CHKERRQ(ierr); 347661095bbSAlp Dener ierr = PetscFree(auglag->Yscatter);CHKERRQ(ierr); 348661095bbSAlp Dener } 349661095bbSAlp Dener } 350661095bbSAlp Dener if (tao->bounded) { 351661095bbSAlp Dener ierr = VecDestroy(&auglag->PL);CHKERRQ(ierr); /* lower bounds for subsolver */ 352661095bbSAlp Dener ierr = VecDestroy(&auglag->PU);CHKERRQ(ierr); /* upper bounds for subsolver */ 353661095bbSAlp Dener } 354661095bbSAlp Dener } 355661095bbSAlp Dener ierr = PetscFree(tao->data);CHKERRQ(ierr); 356661095bbSAlp Dener PetscFunctionReturn(0); 357661095bbSAlp Dener } 358661095bbSAlp Dener 359661095bbSAlp Dener static PetscErrorCode TaoSetFromOptions_ALMM(PetscOptionItems *PetscOptionsObject,Tao tao) 360661095bbSAlp Dener { 361661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM*)tao->data; 362661095bbSAlp Dener PetscInt i; 363661095bbSAlp Dener PetscErrorCode ierr; 364661095bbSAlp Dener 365661095bbSAlp Dener PetscFunctionBegin; 366661095bbSAlp Dener ierr = PetscOptionsHead(PetscOptionsObject,"Augmented Lagrangian multipler method solves problems with general constraints by converting them into a sequence of unconstrained problems.");CHKERRQ(ierr); 367661095bbSAlp Dener ierr = PetscOptionsReal("-tao_almm_mu_init","initial penalty parameter","",auglag->mu0,&auglag->mu0,NULL);CHKERRQ(ierr); 368661095bbSAlp Dener ierr = PetscOptionsReal("-tao_almm_mu_factor","increase factor for the penalty parameter","",auglag->mu_fac,&auglag->mu_fac,NULL);CHKERRQ(ierr); 369661095bbSAlp Dener ierr = PetscOptionsReal("-tao_almm_mu_power_good","exponential for penalty parameter when multiplier update is accepted","",auglag->mu_pow_good,&auglag->mu_pow_good,NULL);CHKERRQ(ierr); 370661095bbSAlp Dener ierr = PetscOptionsReal("-tao_almm_mu_power_bad","exponential for penalty parameter when multiplier update is rejected","",auglag->mu_pow_bad,&auglag->mu_pow_bad,NULL);CHKERRQ(ierr); 371661095bbSAlp Dener ierr = PetscOptionsReal("-tao_almm_mu_max","maximum safeguard for penalty parameter updates","",auglag->mu_max,&auglag->mu_max,NULL);CHKERRQ(ierr); 372661095bbSAlp Dener ierr = PetscOptionsReal("-tao_almm_ye_min","minimum safeguard for equality multiplier updates","",auglag->ye_min,&auglag->ye_min,NULL);CHKERRQ(ierr); 373661095bbSAlp Dener ierr = PetscOptionsReal("-tao_almm_ye_max","maximum safeguard for equality multipliers updates","",auglag->ye_max,&auglag->ye_max,NULL);CHKERRQ(ierr); 374661095bbSAlp Dener ierr = PetscOptionsReal("-tao_almm_yi_min","minimum safeguard for inequality multipliers updates","",auglag->yi_min,&auglag->yi_min,NULL);CHKERRQ(ierr); 375661095bbSAlp Dener ierr = PetscOptionsReal("-tao_almm_yi_max","maximum safeguard for inequality multipliers updates","",auglag->yi_max,&auglag->yi_max,NULL);CHKERRQ(ierr); 376661095bbSAlp Dener ierr = PetscOptionsEnum("-tao_almm_type","augmented Lagrangian formulation type for the subproblem","TaoALMMType",TaoALMMTypes,(PetscEnum)auglag->type,(PetscEnum*)&auglag->type,NULL);CHKERRQ(ierr); 377661095bbSAlp Dener ierr = PetscOptionsTail();CHKERRQ(ierr); 3787721a36fSStefano Zampini ierr = TaoSetOptionsPrefix(auglag->subsolver,((PetscObject)tao)->prefix);CHKERRQ(ierr); 3797721a36fSStefano Zampini ierr = TaoAppendOptionsPrefix(auglag->subsolver,"tao_almm_subsolver_");CHKERRQ(ierr); 380661095bbSAlp Dener ierr = TaoSetFromOptions(auglag->subsolver);CHKERRQ(ierr); 381661095bbSAlp Dener for (i=0; i<tao->numbermonitors; i++) { 382661095bbSAlp Dener ierr = PetscObjectReference((PetscObject)tao->monitorcontext[i]);CHKERRQ(ierr); 383661095bbSAlp Dener ierr = TaoSetMonitor(auglag->subsolver, tao->monitor[i], tao->monitorcontext[i], tao->monitordestroy[i]);CHKERRQ(ierr); 384661095bbSAlp Dener if (tao->monitor[i] == TaoMonitorDefault || tao->monitor[i] == TaoDefaultCMonitor || tao->monitor[i] == TaoDefaultGMonitor || tao->monitor[i] == TaoDefaultSMonitor) { 385661095bbSAlp Dener auglag->info = PETSC_TRUE; 386661095bbSAlp Dener } 387661095bbSAlp Dener } 388661095bbSAlp Dener PetscFunctionReturn(0); 389661095bbSAlp Dener } 390661095bbSAlp Dener 391661095bbSAlp Dener /* -------------------------------------------------------- */ 392661095bbSAlp Dener 393661095bbSAlp Dener /*MC 394661095bbSAlp Dener TaoALMM - Augmented Lagrangian multiplier method for solving nonlinear optimization problems with general constraints. 395661095bbSAlp Dener 396661095bbSAlp Dener Options Database Keys: 397661095bbSAlp Dener + -tao_almm_mu_init <real> - initial penalty parameter (default: 10.) 398661095bbSAlp Dener . -tao_almm_mu_factor <real> - increase factor for the penalty parameter (default: 100.) 399661095bbSAlp Dener . -tao_almm_mu_max <real> - maximum safeguard for penalty parameter updates (default: 1.e20) 400661095bbSAlp Dener . -tao_almm_mu_power_good <real> - exponential for penalty parameter when multiplier update is accepted (default: 0.9) 401661095bbSAlp Dener . -tao_almm_mu_power_bad <real> - exponential for penalty parameter when multiplier update is rejected (default: 0.1) 402661095bbSAlp Dener . -tao_almm_ye_min <real> - minimum safeguard for equality multiplier updates (default: -1.e20) 403661095bbSAlp Dener . -tao_almm_ye_max <real> - maximum safeguard for equality multiplier updates (default: 1.e20) 404661095bbSAlp Dener . -tao_almm_yi_min <real> - minimum safeguard for inequality multiplier updates (default: -1.e20) 405661095bbSAlp Dener . -tao_almm_yi_max <real> - maximum safeguard for inequality multiplier updates (default: 1.e20) 406661095bbSAlp Dener - -tao_almm_type <classic,phr> - change formulation of the augmented Lagrangian merit function for the subproblem (default: classic) 407661095bbSAlp Dener 408661095bbSAlp Dener Level: beginner 409661095bbSAlp Dener 410661095bbSAlp Dener Notes: 411661095bbSAlp Dener This method converts a constrained problem into a sequence of unconstrained problems via the augmented 412661095bbSAlp Dener Lagrangian merit function. Bound constraints are pushed down to the subproblem without any modifications. 413661095bbSAlp Dener 414661095bbSAlp Dener Two formulations are offered for the subproblem: canonical Hestenes-Powell augmented Lagrangian with slack 415661095bbSAlp Dener variables for inequality constraints, and a slack-less Powell-Hestenes-Rockafellar (PHR) formulation utilizing a 416661095bbSAlp Dener pointwise max() penalty on inequality constraints. The canonical augmented Lagrangian formulation typically 417661095bbSAlp Dener converges faster for most problems. However, PHR may be desirable for problems featuring a large number 418661095bbSAlp Dener of inequality constraints because it avoids inflating the size of the subproblem with slack variables. 419661095bbSAlp Dener 420661095bbSAlp Dener The subproblem is solved using a nested first-order TAO solver. The user can retrieve a pointer to 421661095bbSAlp Dener the subsolver via TaoALMMGetSubsolver() or pass command line arguments to it using the 422661095bbSAlp Dener "-tao_almm_subsolver_" prefix. Currently, TaoALMM does not support second-order methods for the 423661095bbSAlp Dener subproblem. It is also highly recommended that the subsolver chosen by the user utilize a trust-region 424661095bbSAlp Dener strategy for globalization (default: TAOBQNKTR) especially if the outer problem features bound constraints. 425661095bbSAlp Dener 426661095bbSAlp Dener .vb 427661095bbSAlp Dener while unconverged 428661095bbSAlp Dener solve argmin_x L(x) s.t. l <= x <= u 429661095bbSAlp Dener if ||c|| <= y_tol 430661095bbSAlp Dener if ||c|| <= c_tol && ||Lgrad|| <= g_tol: 431661095bbSAlp Dener problem converged, return solution 432661095bbSAlp Dener else 433661095bbSAlp Dener constraints sufficiently improved 434661095bbSAlp Dener update multipliers and tighten tolerances 435661095bbSAlp Dener endif 436661095bbSAlp Dener else 437661095bbSAlp Dener constraints did not improve 438661095bbSAlp Dener update penalty and loosen tolerances 439661095bbSAlp Dener endif 440661095bbSAlp Dener endwhile 441661095bbSAlp Dener .ve 442661095bbSAlp Dener 443661095bbSAlp Dener .seealso: TaoALMMGetType(), TaoALMMSetType(), TaoALMMSetSubsolver(), TaoALMMGetSubsolver(), 444661095bbSAlp Dener TaoALMMGetMultipliers(), TaoALMMSetMultipliers(), TaoALMMGetPrimalIS(), TaoALMMGetDualIS() 445661095bbSAlp Dener M*/ 446661095bbSAlp Dener PETSC_EXTERN PetscErrorCode TaoCreate_ALMM(Tao tao) 447661095bbSAlp Dener { 448661095bbSAlp Dener TAO_ALMM *auglag; 449661095bbSAlp Dener PetscErrorCode ierr; 450661095bbSAlp Dener 451661095bbSAlp Dener PetscFunctionBegin; 452661095bbSAlp Dener ierr = PetscNewLog(tao, &auglag);CHKERRQ(ierr); 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 487661095bbSAlp Dener ierr = TaoCreate(PetscObjectComm((PetscObject)tao),&auglag->subsolver);CHKERRQ(ierr); 488661095bbSAlp Dener ierr = TaoSetType(auglag->subsolver, TAOBQNKTR);CHKERRQ(ierr); 489661095bbSAlp Dener ierr = TaoSetTolerances(auglag->subsolver, auglag->gtol, 0.0, 0.0);CHKERRQ(ierr); 490661095bbSAlp Dener ierr = TaoSetMaximumIterations(auglag->subsolver, 1000);CHKERRQ(ierr); 491661095bbSAlp Dener ierr = TaoSetMaximumFunctionEvaluations(auglag->subsolver, 10000);CHKERRQ(ierr); 492661095bbSAlp Dener ierr = TaoSetFunctionLowerBound(auglag->subsolver, PETSC_NINFINITY);CHKERRQ(ierr); 493661095bbSAlp Dener ierr = PetscObjectIncrementTabLevel((PetscObject)auglag->subsolver,(PetscObject)tao,1);CHKERRQ(ierr); 494661095bbSAlp Dener 495661095bbSAlp Dener ierr = PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetType_C", TaoALMMGetType_Private);CHKERRQ(ierr); 496661095bbSAlp Dener ierr = PetscObjectComposeFunction((PetscObject)tao, "TaoALMMSetType_C", TaoALMMSetType_Private);CHKERRQ(ierr); 497661095bbSAlp Dener ierr = PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetSubsolver_C", TaoALMMGetSubsolver_Private);CHKERRQ(ierr); 498661095bbSAlp Dener ierr = PetscObjectComposeFunction((PetscObject)tao, "TaoALMMSetSubsolver_C", TaoALMMSetSubsolver_Private);CHKERRQ(ierr); 499661095bbSAlp Dener ierr = PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetMultipliers_C", TaoALMMGetMultipliers_Private);CHKERRQ(ierr); 500661095bbSAlp Dener ierr = PetscObjectComposeFunction((PetscObject)tao, "TaoALMMSetMultipliers_C", TaoALMMSetMultipliers_Private);CHKERRQ(ierr); 501661095bbSAlp Dener ierr = PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetPrimalIS_C", TaoALMMGetPrimalIS_Private);CHKERRQ(ierr); 502661095bbSAlp Dener ierr = PetscObjectComposeFunction((PetscObject)tao, "TaoALMMGetDualIS_C", TaoALMMGetDualIS_Private);CHKERRQ(ierr); 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 PetscErrorCode ierr; 510661095bbSAlp Dener 511661095bbSAlp Dener PetscFunctionBegin; 512661095bbSAlp Dener if (tao->ineq_constrained) { 513661095bbSAlp Dener ierr = VecScatterBegin(auglag->Pscatter[0], X, P, INSERT_VALUES, SCATTER_REVERSE);CHKERRQ(ierr); 514661095bbSAlp Dener ierr = VecScatterEnd(auglag->Pscatter[0], X, P, INSERT_VALUES, SCATTER_REVERSE);CHKERRQ(ierr); 515661095bbSAlp Dener ierr = VecScatterBegin(auglag->Pscatter[1], S, P, INSERT_VALUES, SCATTER_REVERSE);CHKERRQ(ierr); 516661095bbSAlp Dener ierr = VecScatterEnd(auglag->Pscatter[1], S, P, INSERT_VALUES, SCATTER_REVERSE);CHKERRQ(ierr); 517661095bbSAlp Dener } else { 518661095bbSAlp Dener ierr = VecCopy(X, P);CHKERRQ(ierr); 519661095bbSAlp Dener } 520661095bbSAlp Dener PetscFunctionReturn(0); 521661095bbSAlp Dener } 522661095bbSAlp Dener 523661095bbSAlp Dener static PetscErrorCode TaoALMMCombineDual_Private(Tao tao, Vec EQ, Vec IN, Vec Y) 524661095bbSAlp Dener { 525661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM*)tao->data; 526661095bbSAlp Dener PetscErrorCode ierr; 527661095bbSAlp Dener 528661095bbSAlp Dener PetscFunctionBegin; 529661095bbSAlp Dener if (tao->eq_constrained) { 530661095bbSAlp Dener if (tao->ineq_constrained) { 531661095bbSAlp Dener ierr = VecScatterBegin(auglag->Yscatter[0], EQ, Y, INSERT_VALUES, SCATTER_REVERSE);CHKERRQ(ierr); 532661095bbSAlp Dener ierr = VecScatterEnd(auglag->Yscatter[0], EQ, Y, INSERT_VALUES, SCATTER_REVERSE);CHKERRQ(ierr); 533661095bbSAlp Dener ierr = VecScatterBegin(auglag->Yscatter[1], IN, Y, INSERT_VALUES, SCATTER_REVERSE);CHKERRQ(ierr); 534661095bbSAlp Dener ierr = VecScatterEnd(auglag->Yscatter[1], IN, Y, INSERT_VALUES, SCATTER_REVERSE);CHKERRQ(ierr); 535661095bbSAlp Dener } else { 536661095bbSAlp Dener ierr = VecCopy(EQ, Y);CHKERRQ(ierr); 537661095bbSAlp Dener } 538661095bbSAlp Dener } else { 539661095bbSAlp Dener ierr = VecCopy(IN, Y);CHKERRQ(ierr); 540661095bbSAlp Dener } 541661095bbSAlp Dener PetscFunctionReturn(0); 542661095bbSAlp Dener } 543661095bbSAlp Dener 544661095bbSAlp Dener static PetscErrorCode TaoALMMSplitPrimal_Private(Tao tao, Vec P, Vec X, Vec S) 545661095bbSAlp Dener { 546661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM*)tao->data; 547661095bbSAlp Dener PetscErrorCode ierr; 548661095bbSAlp Dener 549661095bbSAlp Dener PetscFunctionBegin; 550661095bbSAlp Dener if (tao->ineq_constrained) { 551661095bbSAlp Dener ierr = VecScatterBegin(auglag->Pscatter[0], P, X, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr); 552661095bbSAlp Dener ierr = VecScatterEnd(auglag->Pscatter[0], P, X, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr); 553661095bbSAlp Dener ierr = VecScatterBegin(auglag->Pscatter[1], P, S, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr); 554661095bbSAlp Dener ierr = VecScatterEnd(auglag->Pscatter[1], P, S, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr); 555661095bbSAlp Dener } else { 556661095bbSAlp Dener ierr = VecCopy(P, X);CHKERRQ(ierr); 557661095bbSAlp Dener } 558661095bbSAlp Dener PetscFunctionReturn(0); 559661095bbSAlp Dener } 560661095bbSAlp Dener 561661095bbSAlp Dener /* this assumes that the latest constraints are stored in Ce and Ci, and also combined in C */ 562661095bbSAlp Dener static PetscErrorCode TaoALMMComputeOptimalityNorms_Private(Tao tao) 563661095bbSAlp Dener { 564661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM*)tao->data; 565661095bbSAlp Dener PetscErrorCode ierr; 566661095bbSAlp Dener 567661095bbSAlp Dener PetscFunctionBegin; 568661095bbSAlp Dener /* if bounded, project the gradient */ 569661095bbSAlp Dener if (tao->bounded) { 5701e1ea65dSPierre Jolivet ierr = VecBoundGradientProjection(auglag->LgradX, auglag->Px, tao->XL, tao->XU, auglag->LgradX);CHKERRQ(ierr); 571661095bbSAlp Dener } 572661095bbSAlp Dener if (auglag->type == TAO_ALMM_PHR) { 573661095bbSAlp Dener ierr = VecNorm(auglag->LgradX, NORM_INFINITY, &auglag->gnorm);CHKERRQ(ierr); 574661095bbSAlp Dener auglag->cenorm = 0.0; 575661095bbSAlp Dener if (tao->eq_constrained) { 576661095bbSAlp Dener ierr = VecNorm(auglag->Ce, NORM_INFINITY, &auglag->cenorm);CHKERRQ(ierr); 577661095bbSAlp Dener } 578661095bbSAlp Dener auglag->cinorm = 0.0; 579661095bbSAlp Dener if (tao->ineq_constrained) { 580661095bbSAlp Dener ierr = VecCopy(auglag->Yi, auglag->Ciwork);CHKERRQ(ierr); 581661095bbSAlp Dener ierr = VecScale(auglag->Ciwork, -1.0/auglag->mu);CHKERRQ(ierr); 582661095bbSAlp Dener ierr = VecPointwiseMax(auglag->Ciwork, auglag->Ci, auglag->Ciwork);CHKERRQ(ierr); 583661095bbSAlp Dener ierr = VecNorm(auglag->Ciwork, NORM_INFINITY, &auglag->cinorm);CHKERRQ(ierr); 584661095bbSAlp Dener } 585661095bbSAlp Dener auglag->cnorm_old = auglag->cnorm; 586661095bbSAlp Dener auglag->cnorm = PetscMax(auglag->cenorm, auglag->cinorm); 587661095bbSAlp Dener auglag->ytol = auglag->ytol0 * auglag->cnorm_old; 588661095bbSAlp Dener } else { 589661095bbSAlp Dener ierr = VecNorm(auglag->LgradX, NORM_2, &auglag->gnorm);CHKERRQ(ierr); 590661095bbSAlp Dener ierr = VecNorm(auglag->C, NORM_2, &auglag->cnorm);CHKERRQ(ierr); 591661095bbSAlp Dener } 592661095bbSAlp Dener PetscFunctionReturn(0); 593661095bbSAlp Dener } 594661095bbSAlp Dener 595661095bbSAlp Dener static PetscErrorCode TaoALMMEvaluateIterate_Private(Tao tao, Vec P) 596661095bbSAlp Dener { 597661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM*)tao->data; 598661095bbSAlp Dener PetscErrorCode ierr; 599661095bbSAlp Dener 600661095bbSAlp Dener PetscFunctionBegin; 601661095bbSAlp Dener /* split solution into primal and slack components */ 602661095bbSAlp Dener ierr = TaoALMMSplitPrimal_Private(tao, auglag->P, auglag->Px, auglag->Ps);CHKERRQ(ierr); 603661095bbSAlp Dener 604661095bbSAlp Dener /* compute f, df/dx and the constraints */ 605661095bbSAlp Dener ierr = TaoComputeObjectiveAndGradient(tao, auglag->Px, &auglag->fval, auglag->LgradX);CHKERRQ(ierr); 606661095bbSAlp Dener if (tao->eq_constrained) { 607661095bbSAlp Dener ierr = TaoComputeEqualityConstraints(tao, auglag->Px, auglag->Ce);CHKERRQ(ierr); 608661095bbSAlp Dener ierr = TaoComputeJacobianEquality(tao, auglag->Px, auglag->Ae, auglag->Ae);CHKERRQ(ierr); 609661095bbSAlp Dener } 610661095bbSAlp Dener if (tao->ineq_constrained) { 611661095bbSAlp Dener ierr = TaoComputeInequalityConstraints(tao, auglag->Px, auglag->Ci);CHKERRQ(ierr); 612661095bbSAlp Dener ierr = TaoComputeJacobianInequality(tao, auglag->Px, auglag->Ai, auglag->Ai);CHKERRQ(ierr); 613661095bbSAlp Dener switch (auglag->type) { 6147721a36fSStefano Zampini case TAO_ALMM_CLASSIC: 615661095bbSAlp Dener /* classic formulation converts inequality to equality constraints via slack variables */ 616661095bbSAlp Dener ierr = VecAXPY(auglag->Ci, -1.0, auglag->Ps);CHKERRQ(ierr); 617661095bbSAlp Dener break; 6187721a36fSStefano Zampini case TAO_ALMM_PHR: 619661095bbSAlp Dener /* PHR is based on Ci <= 0 while TAO defines Ci >= 0 so we hit it with a negative sign */ 620661095bbSAlp Dener ierr = VecScale(auglag->Ci, -1.0);CHKERRQ(ierr); 621661095bbSAlp Dener ierr = MatScale(auglag->Ai, -1.0);CHKERRQ(ierr); 622661095bbSAlp Dener break; 623661095bbSAlp Dener default: 624661095bbSAlp Dener break; 625661095bbSAlp Dener } 626661095bbSAlp Dener } 627661095bbSAlp Dener /* combine constraints into one vector */ 628661095bbSAlp Dener ierr = TaoALMMCombineDual_Private(tao, auglag->Ce, auglag->Ci, auglag->C);CHKERRQ(ierr); 629661095bbSAlp Dener PetscFunctionReturn(0); 630661095bbSAlp Dener } 631661095bbSAlp Dener 632661095bbSAlp Dener /* 633661095bbSAlp Dener Lphr = f + 0.5*mu*[ (Ce + Ye/mu)^T (Ce + Ye/mu) + pmin(0, Ci + Yi/mu)^T pmin(0, Ci + Yi/mu)] 634661095bbSAlp Dener 635661095bbSAlp Dener dLphr/dX = dF/dX + mu*[ (Ce + Ye/mu)^T Ae + pmin(0, Ci + Yi/mu)^T Ai] 636661095bbSAlp Dener 637661095bbSAlp Dener dLphr/dS = 0 638661095bbSAlp Dener */ 639661095bbSAlp Dener static PetscErrorCode TaoALMMComputePHRLagAndGradient_Private(Tao tao) 640661095bbSAlp Dener { 641661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM*)tao->data; 642661095bbSAlp Dener PetscReal eq_norm=0.0, ineq_norm=0.0; 643661095bbSAlp Dener PetscErrorCode ierr; 644661095bbSAlp Dener 645661095bbSAlp Dener PetscFunctionBegin; 646661095bbSAlp Dener ierr = TaoALMMEvaluateIterate_Private(tao, auglag->P);CHKERRQ(ierr); 647661095bbSAlp Dener if (tao->eq_constrained) { 648661095bbSAlp Dener /* Ce_work = mu*(Ce + Ye/mu) */ 649661095bbSAlp Dener ierr = VecWAXPY(auglag->Cework, 1.0/auglag->mu, auglag->Ye, auglag->Ce);CHKERRQ(ierr); 650661095bbSAlp Dener ierr = VecDot(auglag->Cework, auglag->Cework, &eq_norm);CHKERRQ(ierr); /* contribution to scalar Lagrangian */ 651661095bbSAlp Dener ierr = VecScale(auglag->Cework, auglag->mu);CHKERRQ(ierr); 652661095bbSAlp Dener /* dL/dX += mu*(Ce + Ye/mu)^T Ae */ 653661095bbSAlp Dener ierr = MatMultTransposeAdd(auglag->Ae, auglag->Cework, auglag->LgradX, auglag->LgradX);CHKERRQ(ierr); 654661095bbSAlp Dener } 655661095bbSAlp Dener if (tao->ineq_constrained) { 656661095bbSAlp Dener /* Ci_work = mu * pmax(0, Ci + Yi/mu) where pmax() is pointwise max() */ 657661095bbSAlp Dener ierr = VecWAXPY(auglag->Ciwork, 1.0/auglag->mu, auglag->Yi, auglag->Ci);CHKERRQ(ierr); 658661095bbSAlp Dener ierr = VecPointwiseMax(auglag->Ciwork, auglag->Cizero, auglag->Ciwork);CHKERRQ(ierr); 659661095bbSAlp Dener ierr = VecDot(auglag->Ciwork, auglag->Ciwork, &ineq_norm);CHKERRQ(ierr); /* contribution to scalar Lagrangian */ 660661095bbSAlp Dener /* dL/dX += mu * pmax(0, Ci + Yi/mu)^T Ai */ 661661095bbSAlp Dener ierr = VecScale(auglag->Ciwork, auglag->mu);CHKERRQ(ierr); 662661095bbSAlp Dener ierr = MatMultTransposeAdd(auglag->Ai, auglag->Ciwork, auglag->LgradX, auglag->LgradX);CHKERRQ(ierr); 663a5b23f4aSJose E. Roman /* dL/dS = 0 because there are no slacks in PHR */ 664661095bbSAlp Dener ierr = VecZeroEntries(auglag->LgradS);CHKERRQ(ierr); 665661095bbSAlp Dener } 666661095bbSAlp Dener /* combine gradient together */ 667661095bbSAlp Dener ierr = TaoALMMCombinePrimal_Private(tao, auglag->LgradX, auglag->LgradS, auglag->G);CHKERRQ(ierr); 668661095bbSAlp 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)] */ 669661095bbSAlp Dener auglag->Lval = auglag->fval + 0.5*auglag->mu*(eq_norm + ineq_norm);CHKERRQ(ierr); 670661095bbSAlp Dener PetscFunctionReturn(0); 671661095bbSAlp Dener } 672661095bbSAlp Dener 673661095bbSAlp Dener /* 674661095bbSAlp Dener Lc = F + Ye^TCe + Yi^T(Ci - S) + 0.5*mu*[Ce^TCe + (Ci - S)^T(Ci - S)] 675661095bbSAlp Dener 676661095bbSAlp Dener dLc/dX = dF/dX + Ye^TAe + Yi^TAi + 0.5*mu*[Ce^TAe + (Ci - S)^TAi] 677661095bbSAlp Dener 678661095bbSAlp Dener dLc/dS = -[Yi + mu*(Ci - S)] 679661095bbSAlp Dener */ 680661095bbSAlp Dener static PetscErrorCode TaoALMMComputeAugLagAndGradient_Private(Tao tao) 681661095bbSAlp Dener { 682661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM*)tao->data; 683661095bbSAlp Dener PetscReal yeTce=0.0, yiTcims=0.0, ceTce=0.0, cimsTcims=0.0; 684661095bbSAlp Dener PetscErrorCode ierr; 685661095bbSAlp Dener 686661095bbSAlp Dener PetscFunctionBegin; 687661095bbSAlp Dener ierr = TaoALMMEvaluateIterate_Private(tao, auglag->P);CHKERRQ(ierr); 688661095bbSAlp Dener if (tao->eq_constrained) { 689661095bbSAlp Dener /* compute scalar contributions */ 690661095bbSAlp Dener ierr = VecDot(auglag->Ye, auglag->Ce, &yeTce);CHKERRQ(ierr); 691661095bbSAlp Dener ierr = VecDot(auglag->Ce, auglag->Ce, &ceTce);CHKERRQ(ierr); 692661095bbSAlp Dener /* dL/dX += ye^T Ae */ 693661095bbSAlp Dener ierr = MatMultTransposeAdd(auglag->Ae, auglag->Ye, auglag->LgradX, auglag->LgradX);CHKERRQ(ierr); 694661095bbSAlp Dener /* dL/dX += mu * ce^T Ae */ 695661095bbSAlp Dener ierr = MatMultTranspose(auglag->Ae, auglag->Ce, auglag->Xwork);CHKERRQ(ierr); 696661095bbSAlp Dener ierr = VecAXPY(auglag->LgradX, auglag->mu, auglag->Xwork);CHKERRQ(ierr); 697661095bbSAlp Dener } 698661095bbSAlp Dener if (tao->ineq_constrained) { 699661095bbSAlp Dener /* compute scalar contributions */ 700661095bbSAlp Dener ierr = VecDot(auglag->Yi, auglag->Ci, &yiTcims);CHKERRQ(ierr); 701661095bbSAlp Dener ierr = VecDot(auglag->Ci, auglag->Ci, &cimsTcims);CHKERRQ(ierr); 702661095bbSAlp Dener /* dL/dX += yi^T Ai */ 703661095bbSAlp Dener ierr = MatMultTransposeAdd(auglag->Ai, auglag->Yi, auglag->LgradX, auglag->LgradX);CHKERRQ(ierr); 704661095bbSAlp Dener /* dL/dX += mu * (ci - s)^T Ai */ 705661095bbSAlp Dener ierr = MatMultTranspose(auglag->Ai, auglag->Ci, auglag->Xwork);CHKERRQ(ierr); 706661095bbSAlp Dener ierr = VecAXPY(auglag->LgradX, auglag->mu, auglag->Xwork);CHKERRQ(ierr); 707661095bbSAlp Dener /* dL/dS = -[yi + mu*(ci - s)] */ 708661095bbSAlp Dener ierr = VecWAXPY(auglag->LgradS, auglag->mu, auglag->Ci, auglag->Yi);CHKERRQ(ierr); 709661095bbSAlp Dener ierr = VecScale(auglag->LgradS, -1.0);CHKERRQ(ierr); 710661095bbSAlp Dener } 711661095bbSAlp Dener /* combine gradient together */ 712661095bbSAlp Dener ierr = TaoALMMCombinePrimal_Private(tao, auglag->LgradX, auglag->LgradS, auglag->G);CHKERRQ(ierr); 713661095bbSAlp Dener /* compute L = f + ye^T ce + yi^T (ci - s) + 0.5*mu*||ce||^2 + 0.5*mu*||ci - s||^2 */ 714661095bbSAlp Dener auglag->Lval = auglag->fval + yeTce + yiTcims + 0.5*auglag->mu*(ceTce + cimsTcims); 715661095bbSAlp Dener PetscFunctionReturn(0); 716661095bbSAlp Dener } 717661095bbSAlp Dener 7187721a36fSStefano Zampini PetscErrorCode TaoALMMSubsolverObjective_Private(Tao tao, Vec P, PetscReal *Lval, void *ctx) 7197721a36fSStefano Zampini { 7207721a36fSStefano Zampini TAO_ALMM *auglag = (TAO_ALMM*)ctx; 7217721a36fSStefano Zampini PetscErrorCode ierr; 7227721a36fSStefano Zampini 7237721a36fSStefano Zampini PetscFunctionBegin; 7247721a36fSStefano Zampini ierr = VecCopy(P, auglag->P);CHKERRQ(ierr); 7257721a36fSStefano Zampini ierr = (*auglag->sub_obj)(auglag->parent);CHKERRQ(ierr); 7267721a36fSStefano Zampini *Lval = auglag->Lval; 7277721a36fSStefano Zampini PetscFunctionReturn(0); 7287721a36fSStefano Zampini } 7297721a36fSStefano Zampini 730661095bbSAlp Dener PetscErrorCode TaoALMMSubsolverObjectiveAndGradient_Private(Tao tao, Vec P, PetscReal *Lval, Vec G, void *ctx) 731661095bbSAlp Dener { 732661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM*)ctx; 733661095bbSAlp Dener PetscErrorCode ierr; 734661095bbSAlp Dener 735661095bbSAlp Dener PetscFunctionBegin; 736661095bbSAlp Dener ierr = VecCopy(P, auglag->P);CHKERRQ(ierr); 737661095bbSAlp Dener ierr = (*auglag->sub_obj)(auglag->parent);CHKERRQ(ierr); 738661095bbSAlp Dener ierr = VecCopy(auglag->G, G);CHKERRQ(ierr); 739661095bbSAlp Dener *Lval = auglag->Lval; 740661095bbSAlp Dener PetscFunctionReturn(0); 741661095bbSAlp Dener } 742