1*954e39ddSJose E. Roman #include <../src/tao/constrained/impls/almm/almm.h> /*I "petsctao.h" I*/ 2661095bbSAlp Dener #include <petsctao.h> 3661095bbSAlp Dener #include <petsc/private/petscimpl.h> 4661095bbSAlp Dener #include <petsc/private/vecimpl.h> 5661095bbSAlp Dener 6661095bbSAlp Dener /*@ 76aad120cSJose E. Roman TaoALMMGetType - Retrieve the augmented Lagrangian formulation type for the subproblem. 8661095bbSAlp Dener 9661095bbSAlp Dener Input Parameters: 10c0298470SBarry Smith . tao - the `Tao` context for the `TAOALMM` solver 11661095bbSAlp Dener 12661095bbSAlp Dener Output Parameters: 13661095bbSAlp Dener . type - augmented Lagragrangian type 14661095bbSAlp Dener 15661095bbSAlp Dener Level: advanced 16661095bbSAlp Dener 17c0298470SBarry Smith .seealso: `Tao`, `TAOALMM`, `TaoALMMSetType()`, `TaoALMMType` 18661095bbSAlp Dener @*/ 19d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoALMMGetType(Tao tao, TaoALMMType *type) 20d71ae5a4SJacob Faibussowitsch { 21661095bbSAlp Dener PetscFunctionBegin; 22661095bbSAlp Dener PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 23661095bbSAlp Dener PetscValidPointer(type, 2); 24cac4c232SBarry Smith PetscUseMethod(tao, "TaoALMMGetType_C", (Tao, TaoALMMType *), (tao, type)); 25661095bbSAlp Dener PetscFunctionReturn(0); 26661095bbSAlp Dener } 27661095bbSAlp Dener 28d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoALMMGetType_Private(Tao tao, TaoALMMType *type) 29d71ae5a4SJacob Faibussowitsch { 30661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM *)tao->data; 31661095bbSAlp Dener 32661095bbSAlp Dener PetscFunctionBegin; 33661095bbSAlp Dener *type = auglag->type; 34661095bbSAlp Dener PetscFunctionReturn(0); 35661095bbSAlp Dener } 36661095bbSAlp Dener 37661095bbSAlp Dener /*@ 38661095bbSAlp Dener TaoALMMSetType - Determine the augmented Lagrangian formulation type for the subproblem. 39661095bbSAlp Dener 40661095bbSAlp Dener Input Parameters: 41c0298470SBarry Smith + tao - the Tao context for the `TAOALMM` solver 42661095bbSAlp Dener - type - augmented Lagragrangian type 43661095bbSAlp Dener 44661095bbSAlp Dener Level: advanced 45661095bbSAlp Dener 46c0298470SBarry Smith .seealso: `Tao`, `TAOALMM`, `TaoALMMGetType()`, `TaoALMMType` 47661095bbSAlp Dener @*/ 48d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoALMMSetType(Tao tao, TaoALMMType type) 49d71ae5a4SJacob Faibussowitsch { 50661095bbSAlp Dener PetscFunctionBegin; 51661095bbSAlp Dener PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 52cac4c232SBarry Smith PetscTryMethod(tao, "TaoALMMSetType_C", (Tao, TaoALMMType), (tao, type)); 53661095bbSAlp Dener PetscFunctionReturn(0); 54661095bbSAlp Dener } 55661095bbSAlp Dener 56d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoALMMSetType_Private(Tao tao, TaoALMMType type) 57d71ae5a4SJacob Faibussowitsch { 58661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM *)tao->data; 59661095bbSAlp Dener 60661095bbSAlp Dener PetscFunctionBegin; 613c859ba3SBarry Smith PetscCheck(!tao->setupcalled, PetscObjectComm((PetscObject)tao), PETSC_ERR_ORDER, "TaoALMMSetType() must be called before TaoSetUp()"); 62661095bbSAlp Dener auglag->type = type; 63661095bbSAlp Dener PetscFunctionReturn(0); 64661095bbSAlp Dener } 65661095bbSAlp Dener 66661095bbSAlp Dener /*@ 67c0298470SBarry Smith TaoALMMGetSubsolver - Retrieve the subsolver being used by `TAOALMM`. 68661095bbSAlp Dener 69661095bbSAlp Dener Input Parameters: 70c0298470SBarry Smith . tao - the `Tao` context for the `TAOALMM` solver 71661095bbSAlp Dener 72661095bbSAlp Dener Output Parameter: 73c0298470SBarry Smith . subsolver - the `Tao` context for the subsolver 74661095bbSAlp Dener 75661095bbSAlp Dener Level: advanced 76661095bbSAlp Dener 77c0298470SBarry Smith .seealso: `Tao`, `TAOALMM`, `TaoALMMSetSubsolver()` 78661095bbSAlp Dener @*/ 79d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoALMMGetSubsolver(Tao tao, Tao *subsolver) 80d71ae5a4SJacob Faibussowitsch { 81661095bbSAlp Dener PetscFunctionBegin; 82661095bbSAlp Dener PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 83661095bbSAlp Dener PetscValidPointer(subsolver, 2); 84cac4c232SBarry Smith PetscUseMethod(tao, "TaoALMMGetSubsolver_C", (Tao, Tao *), (tao, subsolver)); 85661095bbSAlp Dener PetscFunctionReturn(0); 86661095bbSAlp Dener } 87661095bbSAlp Dener 88d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoALMMGetSubsolver_Private(Tao tao, Tao *subsolver) 89d71ae5a4SJacob Faibussowitsch { 90661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM *)tao->data; 91661095bbSAlp Dener 92661095bbSAlp Dener PetscFunctionBegin; 93661095bbSAlp Dener *subsolver = auglag->subsolver; 94661095bbSAlp Dener PetscFunctionReturn(0); 95661095bbSAlp Dener } 96661095bbSAlp Dener 97661095bbSAlp Dener /*@ 98c0298470SBarry Smith TaoALMMSetSubsolver - Changes the subsolver inside `TAOALMM` with the user provided one. 99661095bbSAlp Dener 100661095bbSAlp Dener Input Parameters: 101c0298470SBarry Smith + tao - the `Tao` context for the `TAOALMM` solver 102661095bbSAlp Dener - subsolver - the Tao context for the subsolver 103661095bbSAlp Dener 104661095bbSAlp Dener Level: advanced 105661095bbSAlp Dener 106c0298470SBarry Smith Note: 107c0298470SBarry Smith This is not recommended, instead call `TaoALMMGetSubsolver()` and set the type as desired. 108c0298470SBarry Smith 109c0298470SBarry Smith .seealso: `Tao`, `TAOALMM`, `TaoALMMGetSubsolver()` 110661095bbSAlp Dener @*/ 111d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoALMMSetSubsolver(Tao tao, Tao subsolver) 112d71ae5a4SJacob Faibussowitsch { 113661095bbSAlp Dener PetscFunctionBegin; 114661095bbSAlp Dener PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 115661095bbSAlp Dener PetscValidHeaderSpecific(subsolver, TAO_CLASSID, 2); 116cac4c232SBarry Smith PetscTryMethod(tao, "TaoALMMSetSubsolver_C", (Tao, Tao), (tao, subsolver)); 117661095bbSAlp Dener PetscFunctionReturn(0); 118661095bbSAlp Dener } 119661095bbSAlp Dener 120d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoALMMSetSubsolver_Private(Tao tao, Tao subsolver) 121d71ae5a4SJacob Faibussowitsch { 122661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM *)tao->data; 123661095bbSAlp Dener PetscBool compatible; 124661095bbSAlp Dener 125661095bbSAlp Dener PetscFunctionBegin; 126661095bbSAlp Dener if (subsolver == auglag->subsolver) PetscFunctionReturn(0); 127661095bbSAlp Dener if (tao->bounded) { 1289566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompareAny((PetscObject)subsolver, &compatible, TAOSHELL, TAOBNCG, TAOBQNLS, TAOBQNKLS, TAOBQNKTR, TAOBQNKTL, "")); 1293c859ba3SBarry Smith PetscCheck(compatible, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "Subsolver must be a bound-constrained first-order method"); 130661095bbSAlp Dener } else { 1319566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompareAny((PetscObject)subsolver, &compatible, TAOSHELL, TAOCG, TAOLMVM, TAOBNCG, TAOBQNLS, TAOBQNKLS, TAOBQNKTR, TAOBQNKTL, "")); 1323c859ba3SBarry Smith PetscCheck(compatible, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "Subsolver must be a first-order method"); 133661095bbSAlp Dener } 1343c859ba3SBarry Smith PetscCheck(compatible, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "Subsolver must be a first-order method"); 1359566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)subsolver)); 1369566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&auglag->subsolver)); 137661095bbSAlp Dener auglag->subsolver = subsolver; 138661095bbSAlp Dener if (tao->setupcalled) { 1399566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(auglag->subsolver, auglag->P)); 1409566063dSJacob Faibussowitsch PetscCall(TaoSetObjective(auglag->subsolver, TaoALMMSubsolverObjective_Private, (void *)auglag)); 1419566063dSJacob Faibussowitsch PetscCall(TaoSetObjectiveAndGradient(auglag->subsolver, NULL, TaoALMMSubsolverObjectiveAndGradient_Private, (void *)auglag)); 1429566063dSJacob Faibussowitsch PetscCall(TaoSetVariableBounds(auglag->subsolver, auglag->PL, auglag->PU)); 143661095bbSAlp Dener } 144661095bbSAlp Dener PetscFunctionReturn(0); 145661095bbSAlp Dener } 146661095bbSAlp Dener 147661095bbSAlp Dener /*@ 1486aad120cSJose E. Roman TaoALMMGetMultipliers - Retrieve a pointer to the Lagrange multipliers. 149661095bbSAlp Dener 150661095bbSAlp Dener Input Parameters: 151c0298470SBarry Smith . tao - the `Tao` context for the `TAOALMM` solver 152661095bbSAlp Dener 153661095bbSAlp Dener Output Parameters: 154661095bbSAlp Dener . Y - vector of Lagrange multipliers 155661095bbSAlp Dener 156661095bbSAlp Dener Level: advanced 157661095bbSAlp Dener 158661095bbSAlp Dener Notes: 159661095bbSAlp Dener For problems with both equality and inequality constraints, 160661095bbSAlp Dener the multipliers are combined together as Y = (Ye, Yi). Users 161661095bbSAlp Dener can recover copies of the subcomponents using index sets 162c0298470SBarry Smith provided by `TaoALMMGetDualIS()` and use `VecGetSubVector()`. 163661095bbSAlp Dener 164c0298470SBarry Smith .seealso: `TAOALMM`, `Tao`, `TaoALMMSetMultipliers()`, `TaoALMMGetDualIS()` 165661095bbSAlp Dener @*/ 166d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoALMMGetMultipliers(Tao tao, Vec *Y) 167d71ae5a4SJacob Faibussowitsch { 168661095bbSAlp Dener PetscFunctionBegin; 169661095bbSAlp Dener PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 170661095bbSAlp Dener PetscValidPointer(Y, 2); 171cac4c232SBarry Smith PetscUseMethod(tao, "TaoALMMGetMultipliers_C", (Tao, Vec *), (tao, Y)); 172661095bbSAlp Dener PetscFunctionReturn(0); 173661095bbSAlp Dener } 174661095bbSAlp Dener 175d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoALMMGetMultipliers_Private(Tao tao, Vec *Y) 176d71ae5a4SJacob Faibussowitsch { 177661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM *)tao->data; 178661095bbSAlp Dener 179661095bbSAlp Dener PetscFunctionBegin; 1803c859ba3SBarry Smith PetscCheck(tao->setupcalled, PetscObjectComm((PetscObject)tao), PETSC_ERR_ORDER, "TaoSetUp() must be called first for scatters to be constructed"); 181661095bbSAlp Dener *Y = auglag->Y; 182661095bbSAlp Dener PetscFunctionReturn(0); 183661095bbSAlp Dener } 184661095bbSAlp Dener 185661095bbSAlp Dener /*@ 186661095bbSAlp Dener TaoALMMSetMultipliers - Set user-defined Lagrange multipliers. The vector type and 187661095bbSAlp Dener parallel structure of the given vectormust match equality and 188661095bbSAlp Dener inequality constraints. The vector must have a local size equal 189661095bbSAlp Dener to the sum of the local sizes for the constraint vectors, and a 190661095bbSAlp Dener global size equal to the sum of the global sizes of the constraint 191661095bbSAlp Dener vectors. 192661095bbSAlp Dener 193661095bbSAlp Dener Input Parameters: 194c0298470SBarry Smith + tao - the `Tao` context for the `TAOALMM` solver 195661095bbSAlp Dener - Y - vector of Lagrange multipliers 196661095bbSAlp Dener 197661095bbSAlp Dener Level: advanced 198661095bbSAlp Dener 199661095bbSAlp Dener Notes: 200661095bbSAlp Dener This routine is only useful if the user wants to change the 201661095bbSAlp Dener parallel distribution of the combined dual vector in problems that 202661095bbSAlp Dener feature both equality and inequality constraints. For other tasks, 20335cb6cd3SPierre Jolivet it is strongly recommended that the user retrieve the dual vector 204661095bbSAlp Dener created by the solver using TaoALMMGetMultipliers(). 205661095bbSAlp Dener 206c0298470SBarry Smith .seealso: `TAOALMM`, `Tao`, `TaoALMMGetMultipliers()` 207661095bbSAlp Dener @*/ 208d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoALMMSetMultipliers(Tao tao, Vec Y) 209d71ae5a4SJacob Faibussowitsch { 210661095bbSAlp Dener PetscFunctionBegin; 211661095bbSAlp Dener PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 212661095bbSAlp Dener PetscValidHeaderSpecific(Y, VEC_CLASSID, 2); 213cac4c232SBarry Smith PetscTryMethod(tao, "TaoALMMSetMultipliers_C", (Tao, Vec), (tao, Y)); 214661095bbSAlp Dener PetscFunctionReturn(0); 215661095bbSAlp Dener } 216661095bbSAlp Dener 217d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoALMMSetMultipliers_Private(Tao tao, Vec Y) 218d71ae5a4SJacob Faibussowitsch { 219661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM *)tao->data; 220661095bbSAlp Dener VecType Ytype; 221661095bbSAlp Dener PetscInt Nuser, Neq, Nineq, N; 222661095bbSAlp Dener PetscBool same = PETSC_FALSE; 223661095bbSAlp Dener 224661095bbSAlp Dener PetscFunctionBegin; 225661095bbSAlp Dener /* no-op if user provides vector from TaoALMMGetMultipliers() */ 226661095bbSAlp Dener if (Y == auglag->Y) PetscFunctionReturn(0); 227661095bbSAlp Dener /* make sure vector type is same as equality and inequality constraints */ 228661095bbSAlp Dener if (tao->eq_constrained) { 2299566063dSJacob Faibussowitsch PetscCall(VecGetType(tao->constraints_equality, &Ytype)); 230661095bbSAlp Dener } else { 2319566063dSJacob Faibussowitsch PetscCall(VecGetType(tao->constraints_inequality, &Ytype)); 232661095bbSAlp Dener } 2339566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)Y, Ytype, &same)); 2343c859ba3SBarry Smith PetscCheck(same, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "Given vector for multipliers is not the same type as constraint vectors"); 235661095bbSAlp Dener /* make sure global size matches sum of equality and inequality */ 236661095bbSAlp Dener if (tao->eq_constrained) { 2379566063dSJacob Faibussowitsch PetscCall(VecGetSize(tao->constraints_equality, &Neq)); 238661095bbSAlp Dener } else { 239661095bbSAlp Dener Neq = 0; 240661095bbSAlp Dener } 241661095bbSAlp Dener if (tao->ineq_constrained) { 2429566063dSJacob Faibussowitsch PetscCall(VecGetSize(tao->constraints_inequality, &Nineq)); 243661095bbSAlp Dener } else { 244661095bbSAlp Dener Nineq = 0; 245661095bbSAlp Dener } 246661095bbSAlp Dener N = Neq + Nineq; 2479566063dSJacob Faibussowitsch PetscCall(VecGetSize(Y, &Nuser)); 2483c859ba3SBarry Smith PetscCheck(Nuser == N, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "Given vector has wrong global size"); 249661095bbSAlp Dener /* if there is only one type of constraint, then we need the local size to match too */ 250661095bbSAlp Dener if (Neq == 0) { 2519566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(tao->constraints_inequality, &Nineq)); 2529566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(Y, &Nuser)); 2533c859ba3SBarry Smith PetscCheck(Nuser == Nineq, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "Given vector has wrong local size"); 254661095bbSAlp Dener } 255661095bbSAlp Dener if (Nineq == 0) { 2569566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(tao->constraints_equality, &Neq)); 2579566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(Y, &Nuser)); 2583c859ba3SBarry Smith PetscCheck(Nuser == Neq, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "Given vector has wrong local size"); 259661095bbSAlp Dener } 260661095bbSAlp Dener /* if we got here, the given vector is compatible so we can replace the current one */ 2619566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)Y)); 2629566063dSJacob Faibussowitsch PetscCall(VecDestroy(&auglag->Y)); 263661095bbSAlp Dener auglag->Y = Y; 264661095bbSAlp Dener /* if there are both types of constraints and the solver has already been set up, 265661095bbSAlp Dener then we need to regenerate VecScatter objects for the new combined dual vector */ 266661095bbSAlp Dener if (tao->setupcalled && tao->eq_constrained && tao->ineq_constrained) { 2679566063dSJacob Faibussowitsch PetscCall(VecDestroy(&auglag->C)); 2689566063dSJacob Faibussowitsch PetscCall(VecDuplicate(auglag->Y, &auglag->C)); 2699566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&auglag->Yscatter[0])); 2709566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(auglag->Y, auglag->Yis[0], auglag->Ye, NULL, &auglag->Yscatter[0])); 2719566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&auglag->Yscatter[1])); 2729566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(auglag->Y, auglag->Yis[1], auglag->Yi, NULL, &auglag->Yscatter[1])); 273661095bbSAlp Dener } 274661095bbSAlp Dener PetscFunctionReturn(0); 275661095bbSAlp Dener } 276661095bbSAlp Dener 277661095bbSAlp Dener /*@ 2786aad120cSJose E. Roman TaoALMMGetPrimalIS - Retrieve a pointer to the index set that identifies optimization 279661095bbSAlp Dener and slack variable components of the subsolver's solution vector. 280661095bbSAlp Dener Not valid for problems with only equality constraints. 281661095bbSAlp Dener 282f899ff85SJose E. Roman Input Parameter: 283c0298470SBarry Smith . tao - the `Tao` context for the `TAOALMM` solver 284661095bbSAlp Dener 285661095bbSAlp Dener Output Parameters: 286661095bbSAlp Dener + opt_is - index set associated with the optimization variables (NULL if not needed) 287661095bbSAlp Dener - slack_is - index set associated with the slack variables (NULL if not needed) 288661095bbSAlp Dener 289661095bbSAlp Dener Level: advanced 290661095bbSAlp Dener 291c0298470SBarry Smith .seealso: `TAOALMM`, `Tao`, `IS`, `TaoALMMGetPrimalVector()` 292661095bbSAlp Dener @*/ 293d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoALMMGetPrimalIS(Tao tao, IS *opt_is, IS *slack_is) 294d71ae5a4SJacob Faibussowitsch { 295661095bbSAlp Dener PetscFunctionBegin; 296661095bbSAlp Dener PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 297cac4c232SBarry Smith PetscUseMethod(tao, "TaoALMMGetPrimalIS_C", (Tao, IS *, IS *), (tao, opt_is, slack_is)); 298661095bbSAlp Dener PetscFunctionReturn(0); 299661095bbSAlp Dener } 300661095bbSAlp Dener 301d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoALMMGetPrimalIS_Private(Tao tao, IS *opt_is, IS *slack_is) 302d71ae5a4SJacob Faibussowitsch { 303661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM *)tao->data; 304661095bbSAlp Dener 305661095bbSAlp Dener PetscFunctionBegin; 3063c859ba3SBarry Smith PetscCheck(tao->ineq_constrained, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "Primal space has index sets only for inequality constrained problems"); 3073c859ba3SBarry Smith PetscCheck(tao->setupcalled, PetscObjectComm((PetscObject)tao), PETSC_ERR_ORDER, "TaoSetUp() must be called first for index sets to be constructed"); 3087721a36fSStefano Zampini if (opt_is) *opt_is = auglag->Pis[0]; 3097721a36fSStefano Zampini if (slack_is) *slack_is = auglag->Pis[1]; 310661095bbSAlp Dener PetscFunctionReturn(0); 311661095bbSAlp Dener } 312661095bbSAlp Dener 313661095bbSAlp Dener /*@ 3146aad120cSJose E. Roman TaoALMMGetDualIS - Retrieve a pointer to the index set that identifies equality 315661095bbSAlp Dener and inequality constraint components of the dual vector returned 316c0298470SBarry Smith by `TaoALMMGetMultipliers()`. Not valid for problems with only one 317661095bbSAlp Dener type of constraint. 318661095bbSAlp Dener 319f899ff85SJose E. Roman Input Parameter: 320c0298470SBarry Smith . tao - the Tao context for the `TAOALMM` solver 321661095bbSAlp Dener 322661095bbSAlp Dener Output Parameters: 323661095bbSAlp Dener + eq_is - index set associated with the equality constraints (NULL if not needed) 324661095bbSAlp Dener - ineq_is - index set associated with the inequality constraints (NULL if not needed) 325661095bbSAlp Dener 326661095bbSAlp Dener Level: advanced 327661095bbSAlp Dener 328c0298470SBarry Smith .seealso: `TAOALMM`, `Tao`, `TaoALMMGetMultipliers()` 329661095bbSAlp Dener @*/ 330d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoALMMGetDualIS(Tao tao, IS *eq_is, IS *ineq_is) 331d71ae5a4SJacob Faibussowitsch { 332661095bbSAlp Dener PetscFunctionBegin; 333661095bbSAlp Dener PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 334cac4c232SBarry Smith PetscUseMethod(tao, "TaoALMMGetDualIS_C", (Tao, IS *, IS *), (tao, eq_is, ineq_is)); 335661095bbSAlp Dener PetscFunctionReturn(0); 336661095bbSAlp Dener } 337661095bbSAlp Dener 338d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoALMMGetDualIS_Private(Tao tao, IS *eq_is, IS *ineq_is) 339d71ae5a4SJacob Faibussowitsch { 340661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM *)tao->data; 341661095bbSAlp Dener 342661095bbSAlp Dener PetscFunctionBegin; 343661095bbSAlp Dener PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 3443c859ba3SBarry Smith PetscCheck(tao->ineq_constrained && tao->ineq_constrained, PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "Dual space has index sets only when problem has both equality and inequality constraints"); 3453c859ba3SBarry Smith PetscCheck(tao->setupcalled, PetscObjectComm((PetscObject)tao), PETSC_ERR_ORDER, "TaoSetUp() must be called first for index sets to be constructed"); 3467721a36fSStefano Zampini if (eq_is) *eq_is = auglag->Yis[0]; 3477721a36fSStefano Zampini if (ineq_is) *ineq_is = auglag->Yis[1]; 348661095bbSAlp Dener PetscFunctionReturn(0); 349661095bbSAlp Dener } 350