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 /*@ 7661095bbSAlp Dener TaoALMMGetType - Retreive the augmented Lagrangian formulation type for the subproblem. 8661095bbSAlp Dener 9661095bbSAlp Dener Input Parameters: 10661095bbSAlp Dener . 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 17661095bbSAlp Dener .seealso: TAOALMM, TaoALMMSetType(), TaoALMMType 18661095bbSAlp Dener @*/ 19661095bbSAlp Dener PetscErrorCode TaoALMMGetType(Tao tao, TaoALMMType *type) 20661095bbSAlp Dener { 21661095bbSAlp Dener PetscErrorCode ierr; 22661095bbSAlp Dener 23661095bbSAlp Dener PetscFunctionBegin; 24661095bbSAlp Dener PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 25661095bbSAlp Dener PetscValidPointer(type, 2); 26661095bbSAlp Dener ierr = PetscUseMethod(tao,"TaoALMMGetType_C",(Tao,TaoALMMType *),(tao,type));CHKERRQ(ierr); 27661095bbSAlp Dener PetscFunctionReturn(0); 28661095bbSAlp Dener } 29661095bbSAlp Dener 30661095bbSAlp Dener PetscErrorCode TaoALMMGetType_Private(Tao tao, TaoALMMType *type) 31661095bbSAlp Dener { 32661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM*)tao->data; 33661095bbSAlp Dener 34661095bbSAlp Dener PetscFunctionBegin; 35661095bbSAlp Dener *type = auglag->type; 36661095bbSAlp Dener PetscFunctionReturn(0); 37661095bbSAlp Dener } 38661095bbSAlp Dener 39661095bbSAlp Dener /*@ 40661095bbSAlp Dener TaoALMMSetType - Determine the augmented Lagrangian formulation type for the subproblem. 41661095bbSAlp Dener 42661095bbSAlp Dener Input Parameters: 43661095bbSAlp Dener + tao - the Tao context for the TAOALMM solver 44661095bbSAlp Dener - type - augmented Lagragrangian type 45661095bbSAlp Dener 46661095bbSAlp Dener Level: advanced 47661095bbSAlp Dener 48661095bbSAlp Dener .seealso: TAOALMM, TaoALMMGetType(), TaoALMMType 49661095bbSAlp Dener @*/ 50661095bbSAlp Dener PetscErrorCode TaoALMMSetType(Tao tao, TaoALMMType type) 51661095bbSAlp Dener { 52661095bbSAlp Dener PetscErrorCode ierr; 53661095bbSAlp Dener 54661095bbSAlp Dener PetscFunctionBegin; 55661095bbSAlp Dener PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 56661095bbSAlp Dener ierr = PetscTryMethod(tao,"TaoALMMSetType_C",(Tao,TaoALMMType),(tao,type));CHKERRQ(ierr); 57661095bbSAlp Dener PetscFunctionReturn(0); 58661095bbSAlp Dener } 59661095bbSAlp Dener 60661095bbSAlp Dener PetscErrorCode TaoALMMSetType_Private(Tao tao, TaoALMMType type) 61661095bbSAlp Dener { 62661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM*)tao->data; 63661095bbSAlp Dener 64661095bbSAlp Dener PetscFunctionBegin; 65*3c859ba3SBarry Smith PetscCheck(!tao->setupcalled,PetscObjectComm((PetscObject)tao), PETSC_ERR_ORDER, "TaoALMMSetType() must be called before TaoSetUp()"); 66661095bbSAlp Dener auglag->type = type; 67661095bbSAlp Dener PetscFunctionReturn(0); 68661095bbSAlp Dener } 69661095bbSAlp Dener 70661095bbSAlp Dener /*@ 71661095bbSAlp Dener TaoALMMGetSubsolver - Retrieve a pointer to the TAOALMM. 72661095bbSAlp Dener 73661095bbSAlp Dener Input Parameters: 74661095bbSAlp Dener . tao - the Tao context for the TAOALMM solver 75661095bbSAlp Dener 76661095bbSAlp Dener Output Parameter: 77661095bbSAlp Dener . subsolver - the Tao context for the subsolver 78661095bbSAlp Dener 79661095bbSAlp Dener Level: advanced 80661095bbSAlp Dener 81661095bbSAlp Dener .seealso: TAOALMM, TaoALMMSetSubsolver() 82661095bbSAlp Dener @*/ 83661095bbSAlp Dener PetscErrorCode TaoALMMGetSubsolver(Tao tao, Tao *subsolver) 84661095bbSAlp Dener { 85661095bbSAlp Dener PetscErrorCode ierr; 86661095bbSAlp Dener 87661095bbSAlp Dener PetscFunctionBegin; 88661095bbSAlp Dener PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 89661095bbSAlp Dener PetscValidPointer(subsolver, 2); 90661095bbSAlp Dener ierr = PetscUseMethod(tao,"TaoALMMGetSubsolver_C",(Tao,Tao *),(tao,subsolver));CHKERRQ(ierr); 91661095bbSAlp Dener PetscFunctionReturn(0); 92661095bbSAlp Dener } 93661095bbSAlp Dener 94661095bbSAlp Dener PetscErrorCode TaoALMMGetSubsolver_Private(Tao tao, Tao *subsolver) 95661095bbSAlp Dener { 96661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM*)tao->data; 97661095bbSAlp Dener 98661095bbSAlp Dener PetscFunctionBegin; 99661095bbSAlp Dener *subsolver = auglag->subsolver; 100661095bbSAlp Dener PetscFunctionReturn(0); 101661095bbSAlp Dener } 102661095bbSAlp Dener 103661095bbSAlp Dener /*@ 104661095bbSAlp Dener TaoALMMSetSubsolver - Changes the subsolver inside TAOALMM with the user provided one. 105661095bbSAlp Dener 106661095bbSAlp Dener Input Parameters: 107661095bbSAlp Dener + tao - the Tao context for the TAOALMM solver 108661095bbSAlp Dener - subsolver - the Tao context for the subsolver 109661095bbSAlp Dener 110661095bbSAlp Dener Level: advanced 111661095bbSAlp Dener 112661095bbSAlp Dener .seealso: TAOALMM, TaoALMMGetSubsolver() 113661095bbSAlp Dener @*/ 114661095bbSAlp Dener PetscErrorCode TaoALMMSetSubsolver(Tao tao, Tao subsolver) 115661095bbSAlp Dener { 116661095bbSAlp Dener PetscErrorCode ierr; 117661095bbSAlp Dener 118661095bbSAlp Dener PetscFunctionBegin; 119661095bbSAlp Dener PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 120661095bbSAlp Dener PetscValidHeaderSpecific(subsolver, TAO_CLASSID, 2); 121661095bbSAlp Dener ierr = PetscTryMethod(tao,"TaoALMMSetSubsolver_C",(Tao,Tao),(tao,subsolver));CHKERRQ(ierr); 122661095bbSAlp Dener PetscFunctionReturn(0); 123661095bbSAlp Dener } 124661095bbSAlp Dener 125661095bbSAlp Dener PetscErrorCode TaoALMMSetSubsolver_Private(Tao tao, Tao subsolver) 126661095bbSAlp Dener { 127661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM*)tao->data; 128661095bbSAlp Dener PetscBool compatible; 129661095bbSAlp Dener PetscErrorCode ierr; 130661095bbSAlp Dener 131661095bbSAlp Dener PetscFunctionBegin; 132661095bbSAlp Dener if (subsolver == auglag->subsolver) PetscFunctionReturn(0); 133661095bbSAlp Dener if (tao->bounded) { 134661095bbSAlp Dener ierr = PetscObjectTypeCompareAny((PetscObject)subsolver, &compatible, TAOSHELL, TAOBNCG, TAOBQNLS, TAOBQNKLS, TAOBQNKTR, TAOBQNKTL, "");CHKERRQ(ierr); 135*3c859ba3SBarry Smith PetscCheck(compatible,PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "Subsolver must be a bound-constrained first-order method"); 136661095bbSAlp Dener } else { 137661095bbSAlp Dener ierr = PetscObjectTypeCompareAny((PetscObject)subsolver, &compatible, TAOSHELL, TAOCG, TAOLMVM, TAOBNCG, TAOBQNLS, TAOBQNKLS, TAOBQNKTR, TAOBQNKTL, "");CHKERRQ(ierr); 138*3c859ba3SBarry Smith PetscCheck(compatible,PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "Subsolver must be a first-order method"); 139661095bbSAlp Dener } 140*3c859ba3SBarry Smith PetscCheck(compatible,PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "Subsolver must be a first-order method"); 1411e1ea65dSPierre Jolivet ierr = PetscObjectReference((PetscObject)subsolver);CHKERRQ(ierr); 142661095bbSAlp Dener ierr = TaoDestroy(&auglag->subsolver);CHKERRQ(ierr); 143661095bbSAlp Dener auglag->subsolver = subsolver; 144661095bbSAlp Dener if (tao->setupcalled) { 145a82e8c82SStefano Zampini ierr = TaoSetSolution(auglag->subsolver, auglag->P);CHKERRQ(ierr); 146a82e8c82SStefano Zampini ierr = TaoSetObjective(auglag->subsolver, TaoALMMSubsolverObjective_Private, (void*)auglag);CHKERRQ(ierr); 147a82e8c82SStefano Zampini ierr = TaoSetObjectiveAndGradient(auglag->subsolver, NULL, TaoALMMSubsolverObjectiveAndGradient_Private, (void*)auglag);CHKERRQ(ierr); 148661095bbSAlp Dener ierr = TaoSetVariableBounds(auglag->subsolver, auglag->PL, auglag->PU);CHKERRQ(ierr); 149661095bbSAlp Dener } 150661095bbSAlp Dener PetscFunctionReturn(0); 151661095bbSAlp Dener } 152661095bbSAlp Dener 153661095bbSAlp Dener /*@ 154661095bbSAlp Dener TaoALMMGetMultipliers - Retreive a pointer to the Lagrange multipliers. 155661095bbSAlp Dener 156661095bbSAlp Dener Input Parameters: 157661095bbSAlp Dener . tao - the Tao context for the TAOALMM solver 158661095bbSAlp Dener 159661095bbSAlp Dener Output Parameters: 160661095bbSAlp Dener . Y - vector of Lagrange multipliers 161661095bbSAlp Dener 162661095bbSAlp Dener Level: advanced 163661095bbSAlp Dener 164661095bbSAlp Dener Notes: 165661095bbSAlp Dener For problems with both equality and inequality constraints, 166661095bbSAlp Dener the multipliers are combined together as Y = (Ye, Yi). Users 167661095bbSAlp Dener can recover copies of the subcomponents using index sets 168661095bbSAlp Dener provided by TaoALMMGetDualIS() and use VecGetSubVector(). 169661095bbSAlp Dener 170661095bbSAlp Dener .seealso: TAOALMM, TaoALMMSetMultipliers(), TaoALMMGetDualIS() 171661095bbSAlp Dener @*/ 172661095bbSAlp Dener PetscErrorCode TaoALMMGetMultipliers(Tao tao, Vec *Y) 173661095bbSAlp Dener { 174661095bbSAlp Dener PetscErrorCode ierr; 175661095bbSAlp Dener 176661095bbSAlp Dener PetscFunctionBegin; 177661095bbSAlp Dener PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 178661095bbSAlp Dener PetscValidPointer(Y, 2); 179661095bbSAlp Dener ierr = PetscUseMethod(tao,"TaoALMMGetMultipliers_C",(Tao,Vec *),(tao,Y));CHKERRQ(ierr); 180661095bbSAlp Dener PetscFunctionReturn(0); 181661095bbSAlp Dener } 182661095bbSAlp Dener 183661095bbSAlp Dener PetscErrorCode TaoALMMGetMultipliers_Private(Tao tao, Vec *Y) 184661095bbSAlp Dener { 185661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM*)tao->data; 186661095bbSAlp Dener 187661095bbSAlp Dener PetscFunctionBegin; 188*3c859ba3SBarry Smith PetscCheck(tao->setupcalled,PetscObjectComm((PetscObject)tao), PETSC_ERR_ORDER, "TaoSetUp() must be called first for scatters to be constructed"); 189661095bbSAlp Dener *Y = auglag->Y; 190661095bbSAlp Dener PetscFunctionReturn(0); 191661095bbSAlp Dener } 192661095bbSAlp Dener 193661095bbSAlp Dener /*@ 194661095bbSAlp Dener TaoALMMSetMultipliers - Set user-defined Lagrange multipliers. The vector type and 195661095bbSAlp Dener parallel structure of the given vectormust match equality and 196661095bbSAlp Dener inequality constraints. The vector must have a local size equal 197661095bbSAlp Dener to the sum of the local sizes for the constraint vectors, and a 198661095bbSAlp Dener global size equal to the sum of the global sizes of the constraint 199661095bbSAlp Dener vectors. 200661095bbSAlp Dener 201661095bbSAlp Dener Input Parameters: 202661095bbSAlp Dener + tao - the Tao context for the TAOALMM solver 203661095bbSAlp Dener - Y - vector of Lagrange multipliers 204661095bbSAlp Dener 205661095bbSAlp Dener Level: advanced 206661095bbSAlp Dener 207661095bbSAlp Dener Notes: 208661095bbSAlp Dener This routine is only useful if the user wants to change the 209661095bbSAlp Dener parallel distribution of the combined dual vector in problems that 210661095bbSAlp Dener feature both equality and inequality constraints. For other tasks, 211661095bbSAlp Dener it is strongly recommended that the user retreive the dual vector 212661095bbSAlp Dener created by the solver using TaoALMMGetMultipliers(). 213661095bbSAlp Dener 214661095bbSAlp Dener .seealso: TAOALMM, TaoALMMGetMultipliers() 215661095bbSAlp Dener @*/ 216661095bbSAlp Dener PetscErrorCode TaoALMMSetMultipliers(Tao tao, Vec Y) 217661095bbSAlp Dener { 218661095bbSAlp Dener PetscErrorCode ierr; 219661095bbSAlp Dener 220661095bbSAlp Dener PetscFunctionBegin; 221661095bbSAlp Dener PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 222661095bbSAlp Dener PetscValidHeaderSpecific(Y, VEC_CLASSID, 2); 223661095bbSAlp Dener ierr = PetscTryMethod(tao,"TaoALMMSetMultipliers_C",(Tao,Vec),(tao,Y));CHKERRQ(ierr); 224661095bbSAlp Dener PetscFunctionReturn(0); 225661095bbSAlp Dener } 226661095bbSAlp Dener 227661095bbSAlp Dener PetscErrorCode TaoALMMSetMultipliers_Private(Tao tao, Vec Y) 228661095bbSAlp Dener { 229661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM*)tao->data; 230661095bbSAlp Dener VecType Ytype; 231661095bbSAlp Dener PetscInt Nuser, Neq, Nineq, N; 232661095bbSAlp Dener PetscBool same = PETSC_FALSE; 233661095bbSAlp Dener PetscErrorCode ierr; 234661095bbSAlp Dener 235661095bbSAlp Dener PetscFunctionBegin; 236661095bbSAlp Dener /* no-op if user provides vector from TaoALMMGetMultipliers() */ 237661095bbSAlp Dener if (Y == auglag->Y) PetscFunctionReturn(0); 238661095bbSAlp Dener /* make sure vector type is same as equality and inequality constraints */ 239661095bbSAlp Dener if (tao->eq_constrained) { 2401e1ea65dSPierre Jolivet ierr = VecGetType(tao->constraints_equality, &Ytype);CHKERRQ(ierr); 241661095bbSAlp Dener } else { 2421e1ea65dSPierre Jolivet ierr = VecGetType(tao->constraints_inequality, &Ytype);CHKERRQ(ierr); 243661095bbSAlp Dener } 2441e1ea65dSPierre Jolivet ierr = PetscObjectTypeCompare((PetscObject)Y, Ytype, &same);CHKERRQ(ierr); 245*3c859ba3SBarry Smith PetscCheck(same,PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "Given vector for multipliers is not the same type as constraint vectors"); 246661095bbSAlp Dener /* make sure global size matches sum of equality and inequality */ 247661095bbSAlp Dener if (tao->eq_constrained) { 2481e1ea65dSPierre Jolivet ierr = VecGetSize(tao->constraints_equality, &Neq);CHKERRQ(ierr); 249661095bbSAlp Dener } else { 250661095bbSAlp Dener Neq = 0; 251661095bbSAlp Dener } 252661095bbSAlp Dener if (tao->ineq_constrained) { 2531e1ea65dSPierre Jolivet ierr = VecGetSize(tao->constraints_inequality, &Nineq);CHKERRQ(ierr); 254661095bbSAlp Dener } else { 255661095bbSAlp Dener Nineq = 0; 256661095bbSAlp Dener } 257661095bbSAlp Dener N = Neq + Nineq; 258661095bbSAlp Dener ierr = VecGetSize(Y, &Nuser);CHKERRQ(ierr); 259*3c859ba3SBarry Smith PetscCheck(Nuser == N,PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "Given vector has wrong global size"); 260661095bbSAlp Dener /* if there is only one type of constraint, then we need the local size to match too */ 261661095bbSAlp Dener if (Neq == 0) { 262661095bbSAlp Dener ierr = VecGetLocalSize(tao->constraints_inequality, &Nineq);CHKERRQ(ierr); 263661095bbSAlp Dener ierr = VecGetLocalSize(Y, &Nuser);CHKERRQ(ierr); 264*3c859ba3SBarry Smith PetscCheck(Nuser == Nineq,PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "Given vector has wrong local size"); 265661095bbSAlp Dener } 266661095bbSAlp Dener if (Nineq == 0) { 267661095bbSAlp Dener ierr = VecGetLocalSize(tao->constraints_equality, &Neq);CHKERRQ(ierr); 268661095bbSAlp Dener ierr = VecGetLocalSize(Y, &Nuser);CHKERRQ(ierr); 269*3c859ba3SBarry Smith PetscCheck(Nuser == Neq,PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "Given vector has wrong local size"); 270661095bbSAlp Dener } 271661095bbSAlp Dener /* if we got here, the given vector is compatible so we can replace the current one */ 272661095bbSAlp Dener ierr = PetscObjectReference((PetscObject)Y);CHKERRQ(ierr); 273661095bbSAlp Dener ierr = VecDestroy(&auglag->Y);CHKERRQ(ierr); 274661095bbSAlp Dener auglag->Y = Y; 275661095bbSAlp Dener /* if there are both types of constraints and the solver has already been set up, 276661095bbSAlp Dener then we need to regenerate VecScatter objects for the new combined dual vector */ 277661095bbSAlp Dener if (tao->setupcalled && tao->eq_constrained && tao->ineq_constrained) { 278661095bbSAlp Dener ierr = VecDestroy(&auglag->C);CHKERRQ(ierr); 279661095bbSAlp Dener ierr = VecDuplicate(auglag->Y, &auglag->C);CHKERRQ(ierr); 280661095bbSAlp Dener ierr = VecScatterDestroy(&auglag->Yscatter[0]);CHKERRQ(ierr); 281661095bbSAlp Dener ierr = VecScatterCreate(auglag->Y, auglag->Yis[0], auglag->Ye, NULL, &auglag->Yscatter[0]);CHKERRQ(ierr); 282661095bbSAlp Dener ierr = VecScatterDestroy(&auglag->Yscatter[1]);CHKERRQ(ierr); 283661095bbSAlp Dener ierr = VecScatterCreate(auglag->Y, auglag->Yis[1], auglag->Yi, NULL, &auglag->Yscatter[1]);CHKERRQ(ierr); 284661095bbSAlp Dener } 285661095bbSAlp Dener PetscFunctionReturn(0); 286661095bbSAlp Dener } 287661095bbSAlp Dener 288661095bbSAlp Dener /*@ 289661095bbSAlp Dener TaoALMMGetPrimalIS - Retreive a pointer to the index set that identifies optimization 290661095bbSAlp Dener and slack variable components of the subsolver's solution vector. 291661095bbSAlp Dener Not valid for problems with only equality constraints. 292661095bbSAlp Dener 293f899ff85SJose E. Roman Input Parameter: 294661095bbSAlp Dener . tao - the Tao context for the TAOALMM solver 295661095bbSAlp Dener 296661095bbSAlp Dener Output Parameters: 297661095bbSAlp Dener + opt_is - index set associated with the optimization variables (NULL if not needed) 298661095bbSAlp Dener - slack_is - index set associated with the slack variables (NULL if not needed) 299661095bbSAlp Dener 300661095bbSAlp Dener Level: advanced 301661095bbSAlp Dener 302661095bbSAlp Dener .seealso: TAOALMM, TaoALMMGetPrimalVector() 303661095bbSAlp Dener @*/ 304661095bbSAlp Dener PetscErrorCode TaoALMMGetPrimalIS(Tao tao, IS *opt_is, IS *slack_is) 305661095bbSAlp Dener { 306661095bbSAlp Dener PetscErrorCode ierr; 307661095bbSAlp Dener 308661095bbSAlp Dener PetscFunctionBegin; 309661095bbSAlp Dener PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 310661095bbSAlp Dener ierr = PetscUseMethod(tao,"TaoALMMGetPrimalIS_C",(Tao,IS *,IS *),(tao,opt_is,slack_is));CHKERRQ(ierr); 311661095bbSAlp Dener PetscFunctionReturn(0); 312661095bbSAlp Dener } 313661095bbSAlp Dener 314661095bbSAlp Dener PetscErrorCode TaoALMMGetPrimalIS_Private(Tao tao, IS *opt_is, IS *slack_is) 315661095bbSAlp Dener { 316661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM*)tao->data; 317661095bbSAlp Dener 318661095bbSAlp Dener PetscFunctionBegin; 319*3c859ba3SBarry Smith PetscCheck(tao->ineq_constrained,PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "Primal space has index sets only for inequality constrained problems"); 320*3c859ba3SBarry Smith PetscCheck(tao->setupcalled,PetscObjectComm((PetscObject)tao), PETSC_ERR_ORDER, "TaoSetUp() must be called first for index sets to be constructed"); 3217721a36fSStefano Zampini if (opt_is) *opt_is = auglag->Pis[0]; 3227721a36fSStefano Zampini if (slack_is) *slack_is = auglag->Pis[1]; 323661095bbSAlp Dener PetscFunctionReturn(0); 324661095bbSAlp Dener } 325661095bbSAlp Dener 326661095bbSAlp Dener /*@ 327661095bbSAlp Dener TaoALMMGetDualIS - Retreive a pointer to the index set that identifies equality 328661095bbSAlp Dener and inequality constraint components of the dual vector returned 329661095bbSAlp Dener by TaoALMMGetMultipliers(). Not valid for problems with only one 330661095bbSAlp Dener type of constraint. 331661095bbSAlp Dener 332f899ff85SJose E. Roman Input Parameter: 333661095bbSAlp Dener . tao - the Tao context for the TAOALMM solver 334661095bbSAlp Dener 335661095bbSAlp Dener Output Parameters: 336661095bbSAlp Dener + eq_is - index set associated with the equality constraints (NULL if not needed) 337661095bbSAlp Dener - ineq_is - index set associated with the inequality constraints (NULL if not needed) 338661095bbSAlp Dener 339661095bbSAlp Dener Level: advanced 340661095bbSAlp Dener 341661095bbSAlp Dener .seealso: TAOALMM, TaoALMMGetMultipliers() 342661095bbSAlp Dener @*/ 343661095bbSAlp Dener PetscErrorCode TaoALMMGetDualIS(Tao tao, IS *eq_is, IS *ineq_is) 344661095bbSAlp Dener { 345661095bbSAlp Dener PetscErrorCode ierr; 346661095bbSAlp Dener 347661095bbSAlp Dener PetscFunctionBegin; 348661095bbSAlp Dener PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 349661095bbSAlp Dener ierr = PetscUseMethod(tao,"TaoALMMGetDualIS_C",(Tao,IS *,IS *),(tao,eq_is,ineq_is));CHKERRQ(ierr); 350661095bbSAlp Dener PetscFunctionReturn(0); 351661095bbSAlp Dener } 352661095bbSAlp Dener 353661095bbSAlp Dener PetscErrorCode TaoALMMGetDualIS_Private(Tao tao, IS *eq_is, IS *ineq_is) 354661095bbSAlp Dener { 355661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM*)tao->data; 356661095bbSAlp Dener 357661095bbSAlp Dener PetscFunctionBegin; 358661095bbSAlp Dener PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 359*3c859ba3SBarry 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"); 360*3c859ba3SBarry Smith PetscCheck(tao->setupcalled,PetscObjectComm((PetscObject)tao), PETSC_ERR_ORDER, "TaoSetUp() must be called first for index sets to be constructed"); 3617721a36fSStefano Zampini if (eq_is) *eq_is = auglag->Yis[0]; 3627721a36fSStefano Zampini if (ineq_is) *ineq_is = auglag->Yis[1]; 363661095bbSAlp Dener PetscFunctionReturn(0); 364661095bbSAlp Dener } 365