1*661095bbSAlp Dener #include <../src/tao/constrained/impls/almm/almm.h> /*I "petsctao.h" I*/ /*I "petscvec.h" I*/ 2*661095bbSAlp Dener #include <petsctao.h> 3*661095bbSAlp Dener #include <petsc/private/petscimpl.h> 4*661095bbSAlp Dener #include <petsc/private/vecimpl.h> 5*661095bbSAlp Dener 6*661095bbSAlp Dener /*@ 7*661095bbSAlp Dener TaoALMMGetType - Retreive the augmented Lagrangian formulation type for the subproblem. 8*661095bbSAlp Dener 9*661095bbSAlp Dener Input Parameters: 10*661095bbSAlp Dener . tao - the Tao context for the TAOALMM solver 11*661095bbSAlp Dener 12*661095bbSAlp Dener Output Parameters: 13*661095bbSAlp Dener . type - augmented Lagragrangian type 14*661095bbSAlp Dener 15*661095bbSAlp Dener Level: advanced 16*661095bbSAlp Dener 17*661095bbSAlp Dener .seealso: TAOALMM, TaoALMMSetType(), TaoALMMType 18*661095bbSAlp Dener @*/ 19*661095bbSAlp Dener PetscErrorCode TaoALMMGetType(Tao tao, TaoALMMType *type) 20*661095bbSAlp Dener { 21*661095bbSAlp Dener PetscErrorCode ierr; 22*661095bbSAlp Dener 23*661095bbSAlp Dener PetscFunctionBegin; 24*661095bbSAlp Dener PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 25*661095bbSAlp Dener PetscValidPointer(type, 2); 26*661095bbSAlp Dener ierr = PetscUseMethod(tao,"TaoALMMGetType_C",(Tao,TaoALMMType *),(tao,type));CHKERRQ(ierr); 27*661095bbSAlp Dener PetscFunctionReturn(0); 28*661095bbSAlp Dener } 29*661095bbSAlp Dener 30*661095bbSAlp Dener PetscErrorCode TaoALMMGetType_Private(Tao tao, TaoALMMType *type) 31*661095bbSAlp Dener { 32*661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM*)tao->data; 33*661095bbSAlp Dener 34*661095bbSAlp Dener PetscFunctionBegin; 35*661095bbSAlp Dener *type = auglag->type; 36*661095bbSAlp Dener PetscFunctionReturn(0); 37*661095bbSAlp Dener } 38*661095bbSAlp Dener 39*661095bbSAlp Dener /*@ 40*661095bbSAlp Dener TaoALMMSetType - Determine the augmented Lagrangian formulation type for the subproblem. 41*661095bbSAlp Dener 42*661095bbSAlp Dener Input Parameters: 43*661095bbSAlp Dener + tao - the Tao context for the TAOALMM solver 44*661095bbSAlp Dener - type - augmented Lagragrangian type 45*661095bbSAlp Dener 46*661095bbSAlp Dener Level: advanced 47*661095bbSAlp Dener 48*661095bbSAlp Dener .seealso: TAOALMM, TaoALMMGetType(), TaoALMMType 49*661095bbSAlp Dener @*/ 50*661095bbSAlp Dener PetscErrorCode TaoALMMSetType(Tao tao, TaoALMMType type) 51*661095bbSAlp Dener { 52*661095bbSAlp Dener PetscErrorCode ierr; 53*661095bbSAlp Dener 54*661095bbSAlp Dener PetscFunctionBegin; 55*661095bbSAlp Dener PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 56*661095bbSAlp Dener ierr = PetscTryMethod(tao,"TaoALMMSetType_C",(Tao,TaoALMMType),(tao,type));CHKERRQ(ierr); 57*661095bbSAlp Dener PetscFunctionReturn(0); 58*661095bbSAlp Dener } 59*661095bbSAlp Dener 60*661095bbSAlp Dener PetscErrorCode TaoALMMSetType_Private(Tao tao, TaoALMMType type) 61*661095bbSAlp Dener { 62*661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM*)tao->data; 63*661095bbSAlp Dener 64*661095bbSAlp Dener PetscFunctionBegin; 65*661095bbSAlp Dener if (tao->setupcalled) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ORDER, "TaoALMMSetType() must be called before TaoSetUp()"); 66*661095bbSAlp Dener auglag->type = type; 67*661095bbSAlp Dener PetscFunctionReturn(0); 68*661095bbSAlp Dener } 69*661095bbSAlp Dener 70*661095bbSAlp Dener /*@ 71*661095bbSAlp Dener TaoALMMGetSubsolver - Retrieve a pointer to the TAOALMM. 72*661095bbSAlp Dener 73*661095bbSAlp Dener Input Parameters: 74*661095bbSAlp Dener . tao - the Tao context for the TAOALMM solver 75*661095bbSAlp Dener 76*661095bbSAlp Dener Output Parameter: 77*661095bbSAlp Dener . subsolver - the Tao context for the subsolver 78*661095bbSAlp Dener 79*661095bbSAlp Dener Level: advanced 80*661095bbSAlp Dener 81*661095bbSAlp Dener .seealso: TAOALMM, TaoALMMSetSubsolver() 82*661095bbSAlp Dener @*/ 83*661095bbSAlp Dener PetscErrorCode TaoALMMGetSubsolver(Tao tao, Tao *subsolver) 84*661095bbSAlp Dener { 85*661095bbSAlp Dener PetscErrorCode ierr; 86*661095bbSAlp Dener 87*661095bbSAlp Dener PetscFunctionBegin; 88*661095bbSAlp Dener PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 89*661095bbSAlp Dener PetscValidPointer(subsolver, 2); 90*661095bbSAlp Dener ierr = PetscUseMethod(tao,"TaoALMMGetSubsolver_C",(Tao,Tao *),(tao,subsolver));CHKERRQ(ierr); 91*661095bbSAlp Dener PetscFunctionReturn(0); 92*661095bbSAlp Dener } 93*661095bbSAlp Dener 94*661095bbSAlp Dener PetscErrorCode TaoALMMGetSubsolver_Private(Tao tao, Tao *subsolver) 95*661095bbSAlp Dener { 96*661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM*)tao->data; 97*661095bbSAlp Dener 98*661095bbSAlp Dener PetscFunctionBegin; 99*661095bbSAlp Dener *subsolver = auglag->subsolver; 100*661095bbSAlp Dener PetscFunctionReturn(0); 101*661095bbSAlp Dener } 102*661095bbSAlp Dener 103*661095bbSAlp Dener /*@ 104*661095bbSAlp Dener TaoALMMSetSubsolver - Changes the subsolver inside TAOALMM with the user provided one. 105*661095bbSAlp Dener 106*661095bbSAlp Dener Input Parameters: 107*661095bbSAlp Dener + tao - the Tao context for the TAOALMM solver 108*661095bbSAlp Dener - subsolver - the Tao context for the subsolver 109*661095bbSAlp Dener 110*661095bbSAlp Dener Level: advanced 111*661095bbSAlp Dener 112*661095bbSAlp Dener .seealso: TAOALMM, TaoALMMGetSubsolver() 113*661095bbSAlp Dener @*/ 114*661095bbSAlp Dener PetscErrorCode TaoALMMSetSubsolver(Tao tao, Tao subsolver) 115*661095bbSAlp Dener { 116*661095bbSAlp Dener PetscErrorCode ierr; 117*661095bbSAlp Dener 118*661095bbSAlp Dener PetscFunctionBegin; 119*661095bbSAlp Dener PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 120*661095bbSAlp Dener PetscValidHeaderSpecific(subsolver, TAO_CLASSID, 2); 121*661095bbSAlp Dener ierr = PetscTryMethod(tao,"TaoALMMSetSubsolver_C",(Tao,Tao),(tao,subsolver));CHKERRQ(ierr); 122*661095bbSAlp Dener PetscFunctionReturn(0); 123*661095bbSAlp Dener } 124*661095bbSAlp Dener 125*661095bbSAlp Dener PetscErrorCode TaoALMMSetSubsolver_Private(Tao tao, Tao subsolver) 126*661095bbSAlp Dener { 127*661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM*)tao->data; 128*661095bbSAlp Dener PetscBool compatible; 129*661095bbSAlp Dener PetscErrorCode ierr; 130*661095bbSAlp Dener 131*661095bbSAlp Dener PetscFunctionBegin; 132*661095bbSAlp Dener if (subsolver == auglag->subsolver) PetscFunctionReturn(0); 133*661095bbSAlp Dener if (tao->bounded) { 134*661095bbSAlp Dener ierr = PetscObjectTypeCompareAny((PetscObject)subsolver, &compatible, TAOSHELL, TAOBNCG, TAOBQNLS, TAOBQNKLS, TAOBQNKTR, TAOBQNKTL, "");CHKERRQ(ierr); 135*661095bbSAlp Dener if (!compatible) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "Subsolver must be a bound-constrained first-order method"); 136*661095bbSAlp Dener } else { 137*661095bbSAlp Dener ierr = PetscObjectTypeCompareAny((PetscObject)subsolver, &compatible, TAOSHELL, TAOCG, TAOLMVM, TAOBNCG, TAOBQNLS, TAOBQNKLS, TAOBQNKTR, TAOBQNKTL, "");CHKERRQ(ierr); 138*661095bbSAlp Dener if (!compatible) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "Subsolver must be a first-order method"); 139*661095bbSAlp Dener } 140*661095bbSAlp Dener if (!compatible) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "Subsolver must be a first-order method"); 141*661095bbSAlp Dener ierr = PetscObjectReference((PetscObject)subsolver); 142*661095bbSAlp Dener ierr = TaoDestroy(&auglag->subsolver);CHKERRQ(ierr); 143*661095bbSAlp Dener auglag->subsolver = subsolver; 144*661095bbSAlp Dener if (tao->setupcalled) { 145*661095bbSAlp Dener ierr = TaoSetInitialVector(auglag->subsolver, auglag->P);CHKERRQ(ierr); 146*661095bbSAlp Dener ierr = TaoSetObjectiveAndGradientRoutine(auglag->subsolver, TaoALMMSubsolverObjectiveAndGradient_Private, (void*)auglag);CHKERRQ(ierr); 147*661095bbSAlp Dener ierr = TaoSetVariableBounds(auglag->subsolver, auglag->PL, auglag->PU);CHKERRQ(ierr); 148*661095bbSAlp Dener } 149*661095bbSAlp Dener PetscFunctionReturn(0); 150*661095bbSAlp Dener } 151*661095bbSAlp Dener 152*661095bbSAlp Dener /*@ 153*661095bbSAlp Dener TaoALMMGetMultipliers - Retreive a pointer to the Lagrange multipliers. 154*661095bbSAlp Dener 155*661095bbSAlp Dener Input Parameters: 156*661095bbSAlp Dener . tao - the Tao context for the TAOALMM solver 157*661095bbSAlp Dener 158*661095bbSAlp Dener Output Parameters: 159*661095bbSAlp Dener . Y - vector of Lagrange multipliers 160*661095bbSAlp Dener 161*661095bbSAlp Dener Level: advanced 162*661095bbSAlp Dener 163*661095bbSAlp Dener Notes: 164*661095bbSAlp Dener For problems with both equality and inequality constraints, 165*661095bbSAlp Dener the multipliers are combined together as Y = (Ye, Yi). Users 166*661095bbSAlp Dener can recover copies of the subcomponents using index sets 167*661095bbSAlp Dener provided by TaoALMMGetDualIS() and use VecGetSubVector(). 168*661095bbSAlp Dener 169*661095bbSAlp Dener .seealso: TAOALMM, TaoALMMSetMultipliers(), TaoALMMGetDualIS() 170*661095bbSAlp Dener @*/ 171*661095bbSAlp Dener PetscErrorCode TaoALMMGetMultipliers(Tao tao, Vec *Y) 172*661095bbSAlp Dener { 173*661095bbSAlp Dener PetscErrorCode ierr; 174*661095bbSAlp Dener 175*661095bbSAlp Dener PetscFunctionBegin; 176*661095bbSAlp Dener PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 177*661095bbSAlp Dener PetscValidPointer(Y, 2); 178*661095bbSAlp Dener ierr = PetscUseMethod(tao,"TaoALMMGetMultipliers_C",(Tao,Vec *),(tao,Y));CHKERRQ(ierr); 179*661095bbSAlp Dener PetscFunctionReturn(0); 180*661095bbSAlp Dener } 181*661095bbSAlp Dener 182*661095bbSAlp Dener PetscErrorCode TaoALMMGetMultipliers_Private(Tao tao, Vec *Y) 183*661095bbSAlp Dener { 184*661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM*)tao->data; 185*661095bbSAlp Dener 186*661095bbSAlp Dener PetscFunctionBegin; 187*661095bbSAlp Dener if (!tao->setupcalled) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ORDER, "TaoSetUp() must be called first for scatters to be constructed"); 188*661095bbSAlp Dener *Y = auglag->Y; 189*661095bbSAlp Dener PetscFunctionReturn(0); 190*661095bbSAlp Dener } 191*661095bbSAlp Dener 192*661095bbSAlp Dener /*@ 193*661095bbSAlp Dener TaoALMMSetMultipliers - Set user-defined Lagrange multipliers. The vector type and 194*661095bbSAlp Dener parallel structure of the given vectormust match equality and 195*661095bbSAlp Dener inequality constraints. The vector must have a local size equal 196*661095bbSAlp Dener to the sum of the local sizes for the constraint vectors, and a 197*661095bbSAlp Dener global size equal to the sum of the global sizes of the constraint 198*661095bbSAlp Dener vectors. 199*661095bbSAlp Dener 200*661095bbSAlp Dener Input Parameters: 201*661095bbSAlp Dener + tao - the Tao context for the TAOALMM solver 202*661095bbSAlp Dener - Y - vector of Lagrange multipliers 203*661095bbSAlp Dener 204*661095bbSAlp Dener Level: advanced 205*661095bbSAlp Dener 206*661095bbSAlp Dener Notes: 207*661095bbSAlp Dener This routine is only useful if the user wants to change the 208*661095bbSAlp Dener parallel distribution of the combined dual vector in problems that 209*661095bbSAlp Dener feature both equality and inequality constraints. For other tasks, 210*661095bbSAlp Dener it is strongly recommended that the user retreive the dual vector 211*661095bbSAlp Dener created by the solver using TaoALMMGetMultipliers(). 212*661095bbSAlp Dener 213*661095bbSAlp Dener .seealso: TAOALMM, TaoALMMGetMultipliers() 214*661095bbSAlp Dener @*/ 215*661095bbSAlp Dener PetscErrorCode TaoALMMSetMultipliers(Tao tao, Vec Y) 216*661095bbSAlp Dener { 217*661095bbSAlp Dener PetscErrorCode ierr; 218*661095bbSAlp Dener 219*661095bbSAlp Dener PetscFunctionBegin; 220*661095bbSAlp Dener PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 221*661095bbSAlp Dener PetscValidHeaderSpecific(Y, VEC_CLASSID, 2); 222*661095bbSAlp Dener ierr = PetscTryMethod(tao,"TaoALMMSetMultipliers_C",(Tao,Vec),(tao,Y));CHKERRQ(ierr); 223*661095bbSAlp Dener PetscFunctionReturn(0); 224*661095bbSAlp Dener } 225*661095bbSAlp Dener 226*661095bbSAlp Dener PetscErrorCode TaoALMMSetMultipliers_Private(Tao tao, Vec Y) 227*661095bbSAlp Dener { 228*661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM*)tao->data; 229*661095bbSAlp Dener VecType Ytype; 230*661095bbSAlp Dener PetscInt Nuser, Neq, Nineq, N; 231*661095bbSAlp Dener PetscBool same = PETSC_FALSE; 232*661095bbSAlp Dener PetscErrorCode ierr; 233*661095bbSAlp Dener 234*661095bbSAlp Dener PetscFunctionBegin; 235*661095bbSAlp Dener /* no-op if user provides vector from TaoALMMGetMultipliers() */ 236*661095bbSAlp Dener if (Y == auglag->Y) PetscFunctionReturn(0); 237*661095bbSAlp Dener /* make sure vector type is same as equality and inequality constraints */ 238*661095bbSAlp Dener if (tao->eq_constrained) { 239*661095bbSAlp Dener ierr = VecGetType(tao->constraints_equality, &Ytype); 240*661095bbSAlp Dener } else { 241*661095bbSAlp Dener ierr = VecGetType(tao->constraints_inequality, &Ytype); 242*661095bbSAlp Dener } 243*661095bbSAlp Dener ierr = PetscObjectTypeCompare((PetscObject)Y, Ytype, &same); 244*661095bbSAlp Dener if (!same) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "Given vector for multipliers is not the same type as constraint vectors"); 245*661095bbSAlp Dener /* make sure global size matches sum of equality and inequality */ 246*661095bbSAlp Dener if (tao->eq_constrained) { 247*661095bbSAlp Dener ierr = VecGetSize(tao->constraints_equality, &Neq); 248*661095bbSAlp Dener } else { 249*661095bbSAlp Dener Neq = 0; 250*661095bbSAlp Dener } 251*661095bbSAlp Dener if (tao->ineq_constrained) { 252*661095bbSAlp Dener ierr = VecGetSize(tao->constraints_inequality, &Nineq); 253*661095bbSAlp Dener } else { 254*661095bbSAlp Dener Nineq = 0; 255*661095bbSAlp Dener } 256*661095bbSAlp Dener N = Neq + Nineq; 257*661095bbSAlp Dener ierr = VecGetSize(Y, &Nuser);CHKERRQ(ierr); 258*661095bbSAlp Dener if (Nuser != N) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "Given vector has wrong global size"); 259*661095bbSAlp Dener /* if there is only one type of constraint, then we need the local size to match too */ 260*661095bbSAlp Dener if (Neq == 0) { 261*661095bbSAlp Dener ierr = VecGetLocalSize(tao->constraints_inequality, &Nineq);CHKERRQ(ierr); 262*661095bbSAlp Dener ierr = VecGetLocalSize(Y, &Nuser);CHKERRQ(ierr); 263*661095bbSAlp Dener if (Nuser != Nineq) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "Given vector has wrong local size"); 264*661095bbSAlp Dener } 265*661095bbSAlp Dener if (Nineq == 0) { 266*661095bbSAlp Dener ierr = VecGetLocalSize(tao->constraints_equality, &Neq);CHKERRQ(ierr); 267*661095bbSAlp Dener ierr = VecGetLocalSize(Y, &Nuser);CHKERRQ(ierr); 268*661095bbSAlp Dener if (Nuser != Neq) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "Given vector has wrong local size"); 269*661095bbSAlp Dener } 270*661095bbSAlp Dener /* if we got here, the given vector is compatible so we can replace the current one */ 271*661095bbSAlp Dener ierr = PetscObjectReference((PetscObject)Y);CHKERRQ(ierr); 272*661095bbSAlp Dener ierr = VecDestroy(&auglag->Y);CHKERRQ(ierr); 273*661095bbSAlp Dener auglag->Y = Y; 274*661095bbSAlp Dener /* if there are both types of constraints and the solver has already been set up, 275*661095bbSAlp Dener then we need to regenerate VecScatter objects for the new combined dual vector */ 276*661095bbSAlp Dener if (tao->setupcalled && tao->eq_constrained && tao->ineq_constrained) { 277*661095bbSAlp Dener ierr = VecDestroy(&auglag->C);CHKERRQ(ierr); 278*661095bbSAlp Dener ierr = VecDuplicate(auglag->Y, &auglag->C);CHKERRQ(ierr); 279*661095bbSAlp Dener ierr = VecScatterDestroy(&auglag->Yscatter[0]);CHKERRQ(ierr); 280*661095bbSAlp Dener ierr = VecScatterCreate(auglag->Y, auglag->Yis[0], auglag->Ye, NULL, &auglag->Yscatter[0]);CHKERRQ(ierr); 281*661095bbSAlp Dener ierr = VecScatterDestroy(&auglag->Yscatter[1]);CHKERRQ(ierr); 282*661095bbSAlp Dener ierr = VecScatterCreate(auglag->Y, auglag->Yis[1], auglag->Yi, NULL, &auglag->Yscatter[1]);CHKERRQ(ierr); 283*661095bbSAlp Dener } 284*661095bbSAlp Dener PetscFunctionReturn(0); 285*661095bbSAlp Dener } 286*661095bbSAlp Dener 287*661095bbSAlp Dener /*@ 288*661095bbSAlp Dener TaoALMMGetPrimalIS - Retreive a pointer to the index set that identifies optimization 289*661095bbSAlp Dener and slack variable components of the subsolver's solution vector. 290*661095bbSAlp Dener Not valid for problems with only equality constraints. 291*661095bbSAlp Dener 292*661095bbSAlp Dener Input Parameters: 293*661095bbSAlp Dener . tao - the Tao context for the TAOALMM solver 294*661095bbSAlp Dener 295*661095bbSAlp Dener Output Parameters: 296*661095bbSAlp Dener + opt_is - index set associated with the optimization variables (NULL if not needed) 297*661095bbSAlp Dener - slack_is - index set associated with the slack variables (NULL if not needed) 298*661095bbSAlp Dener 299*661095bbSAlp Dener Level: advanced 300*661095bbSAlp Dener 301*661095bbSAlp Dener .seealso: TAOALMM, TaoALMMGetPrimalVector() 302*661095bbSAlp Dener @*/ 303*661095bbSAlp Dener PetscErrorCode TaoALMMGetPrimalIS(Tao tao, IS *opt_is, IS *slack_is) 304*661095bbSAlp Dener { 305*661095bbSAlp Dener PetscErrorCode ierr; 306*661095bbSAlp Dener 307*661095bbSAlp Dener PetscFunctionBegin; 308*661095bbSAlp Dener PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 309*661095bbSAlp Dener ierr = PetscUseMethod(tao,"TaoALMMGetPrimalIS_C",(Tao,IS *,IS *),(tao,opt_is,slack_is));CHKERRQ(ierr); 310*661095bbSAlp Dener PetscFunctionReturn(0); 311*661095bbSAlp Dener } 312*661095bbSAlp Dener 313*661095bbSAlp Dener PetscErrorCode TaoALMMGetPrimalIS_Private(Tao tao, IS *opt_is, IS *slack_is) 314*661095bbSAlp Dener { 315*661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM*)tao->data; 316*661095bbSAlp Dener 317*661095bbSAlp Dener PetscFunctionBegin; 318*661095bbSAlp Dener if (!tao->ineq_constrained) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "Primal space has index sets only for inequality constrained problems"); 319*661095bbSAlp Dener if (!tao->setupcalled) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ORDER, "TaoSetUp() must be called first for index sets to be constructed"); 320*661095bbSAlp Dener if (!opt_is && !slack_is) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_NULL, "Both index set pointers cannot be NULL"); 321*661095bbSAlp Dener if (opt_is) { 322*661095bbSAlp Dener *opt_is = auglag->Pis[0]; 323*661095bbSAlp Dener } 324*661095bbSAlp Dener if (slack_is) { 325*661095bbSAlp Dener *slack_is = auglag->Pis[1]; 326*661095bbSAlp Dener } 327*661095bbSAlp Dener PetscFunctionReturn(0); 328*661095bbSAlp Dener } 329*661095bbSAlp Dener 330*661095bbSAlp Dener /*@ 331*661095bbSAlp Dener TaoALMMGetDualIS - Retreive a pointer to the index set that identifies equality 332*661095bbSAlp Dener and inequality constraint components of the dual vector returned 333*661095bbSAlp Dener by TaoALMMGetMultipliers(). Not valid for problems with only one 334*661095bbSAlp Dener type of constraint. 335*661095bbSAlp Dener 336*661095bbSAlp Dener Input Parameters: 337*661095bbSAlp Dener . tao - the Tao context for the TAOALMM solver 338*661095bbSAlp Dener 339*661095bbSAlp Dener Output Parameters: 340*661095bbSAlp Dener + eq_is - index set associated with the equality constraints (NULL if not needed) 341*661095bbSAlp Dener - ineq_is - index set associated with the inequality constraints (NULL if not needed) 342*661095bbSAlp Dener 343*661095bbSAlp Dener Level: advanced 344*661095bbSAlp Dener 345*661095bbSAlp Dener .seealso: TAOALMM, TaoALMMGetMultipliers() 346*661095bbSAlp Dener @*/ 347*661095bbSAlp Dener PetscErrorCode TaoALMMGetDualIS(Tao tao, IS *eq_is, IS *ineq_is) 348*661095bbSAlp Dener { 349*661095bbSAlp Dener PetscErrorCode ierr; 350*661095bbSAlp Dener 351*661095bbSAlp Dener PetscFunctionBegin; 352*661095bbSAlp Dener PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 353*661095bbSAlp Dener ierr = PetscUseMethod(tao,"TaoALMMGetDualIS_C",(Tao,IS *,IS *),(tao,eq_is,ineq_is));CHKERRQ(ierr); 354*661095bbSAlp Dener PetscFunctionReturn(0); 355*661095bbSAlp Dener } 356*661095bbSAlp Dener 357*661095bbSAlp Dener PetscErrorCode TaoALMMGetDualIS_Private(Tao tao, IS *eq_is, IS *ineq_is) 358*661095bbSAlp Dener { 359*661095bbSAlp Dener TAO_ALMM *auglag = (TAO_ALMM*)tao->data; 360*661095bbSAlp Dener 361*661095bbSAlp Dener PetscFunctionBegin; 362*661095bbSAlp Dener PetscValidHeaderSpecific(tao, TAO_CLASSID, 1); 363*661095bbSAlp Dener if (!tao->ineq_constrained || !tao->ineq_constrained) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "Dual space has index sets only when problem has both equality and inequality constraints"); 364*661095bbSAlp Dener if (!tao->setupcalled) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ORDER, "TaoSetUp() must be called first for index sets to be constructed"); 365*661095bbSAlp Dener if (!eq_is && !ineq_is) SETERRQ(PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_NULL, "Both index set pointers cannot be NULL"); 366*661095bbSAlp Dener if (eq_is) { 367*661095bbSAlp Dener *eq_is = auglag->Yis[0]; 368*661095bbSAlp Dener } 369*661095bbSAlp Dener if (ineq_is) { 370*661095bbSAlp Dener *ineq_is = auglag->Yis[1]; 371*661095bbSAlp Dener } 372*661095bbSAlp Dener PetscFunctionReturn(0); 373*661095bbSAlp Dener }