xref: /petsc/src/tao/constrained/impls/almm/almmutils.c (revision cac4c232dc4f93991e342196e27ef7b0655dac7b)
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   PetscFunctionBegin;
22661095bbSAlp Dener   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
23661095bbSAlp Dener   PetscValidPointer(type, 2);
24*cac4c232SBarry Smith   PetscUseMethod(tao,"TaoALMMGetType_C",(Tao,TaoALMMType *),(tao,type));
25661095bbSAlp Dener   PetscFunctionReturn(0);
26661095bbSAlp Dener }
27661095bbSAlp Dener 
28661095bbSAlp Dener PetscErrorCode TaoALMMGetType_Private(Tao tao, TaoALMMType *type)
29661095bbSAlp Dener {
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:
41661095bbSAlp Dener +  tao - the Tao context for the TAOALMM solver
42661095bbSAlp Dener -  type - augmented Lagragrangian type
43661095bbSAlp Dener 
44661095bbSAlp Dener    Level: advanced
45661095bbSAlp Dener 
46661095bbSAlp Dener .seealso: TAOALMM, TaoALMMGetType(), TaoALMMType
47661095bbSAlp Dener @*/
48661095bbSAlp Dener PetscErrorCode TaoALMMSetType(Tao tao, TaoALMMType type)
49661095bbSAlp Dener {
50661095bbSAlp Dener   PetscFunctionBegin;
51661095bbSAlp Dener   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
52*cac4c232SBarry Smith   PetscTryMethod(tao,"TaoALMMSetType_C",(Tao,TaoALMMType),(tao,type));
53661095bbSAlp Dener   PetscFunctionReturn(0);
54661095bbSAlp Dener }
55661095bbSAlp Dener 
56661095bbSAlp Dener PetscErrorCode TaoALMMSetType_Private(Tao tao, TaoALMMType type)
57661095bbSAlp Dener {
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 /*@
67661095bbSAlp Dener    TaoALMMGetSubsolver - Retrieve a pointer to the TAOALMM.
68661095bbSAlp Dener 
69661095bbSAlp Dener    Input Parameters:
70661095bbSAlp Dener .  tao - the Tao context for the TAOALMM solver
71661095bbSAlp Dener 
72661095bbSAlp Dener    Output Parameter:
73661095bbSAlp Dener .  subsolver - the Tao context for the subsolver
74661095bbSAlp Dener 
75661095bbSAlp Dener    Level: advanced
76661095bbSAlp Dener 
77661095bbSAlp Dener .seealso: TAOALMM, TaoALMMSetSubsolver()
78661095bbSAlp Dener @*/
79661095bbSAlp Dener PetscErrorCode TaoALMMGetSubsolver(Tao tao, Tao *subsolver)
80661095bbSAlp Dener {
81661095bbSAlp Dener   PetscFunctionBegin;
82661095bbSAlp Dener   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
83661095bbSAlp Dener   PetscValidPointer(subsolver, 2);
84*cac4c232SBarry Smith   PetscUseMethod(tao,"TaoALMMGetSubsolver_C",(Tao,Tao *),(tao,subsolver));
85661095bbSAlp Dener   PetscFunctionReturn(0);
86661095bbSAlp Dener }
87661095bbSAlp Dener 
88661095bbSAlp Dener PetscErrorCode TaoALMMGetSubsolver_Private(Tao tao, Tao *subsolver)
89661095bbSAlp Dener {
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 /*@
98661095bbSAlp Dener    TaoALMMSetSubsolver - Changes the subsolver inside TAOALMM with the user provided one.
99661095bbSAlp Dener 
100661095bbSAlp Dener    Input Parameters:
101661095bbSAlp Dener +  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 
106661095bbSAlp Dener .seealso: TAOALMM, TaoALMMGetSubsolver()
107661095bbSAlp Dener @*/
108661095bbSAlp Dener PetscErrorCode TaoALMMSetSubsolver(Tao tao, Tao subsolver)
109661095bbSAlp Dener {
110661095bbSAlp Dener    PetscFunctionBegin;
111661095bbSAlp Dener    PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
112661095bbSAlp Dener    PetscValidHeaderSpecific(subsolver, TAO_CLASSID, 2);
113*cac4c232SBarry Smith    PetscTryMethod(tao,"TaoALMMSetSubsolver_C",(Tao,Tao),(tao,subsolver));
114661095bbSAlp Dener    PetscFunctionReturn(0);
115661095bbSAlp Dener }
116661095bbSAlp Dener 
117661095bbSAlp Dener PetscErrorCode TaoALMMSetSubsolver_Private(Tao tao, Tao subsolver)
118661095bbSAlp Dener {
119661095bbSAlp Dener   TAO_ALMM       *auglag = (TAO_ALMM*)tao->data;
120661095bbSAlp Dener   PetscBool      compatible;
121661095bbSAlp Dener 
122661095bbSAlp Dener   PetscFunctionBegin;
123661095bbSAlp Dener   if (subsolver == auglag->subsolver) PetscFunctionReturn(0);
124661095bbSAlp Dener   if (tao->bounded) {
1259566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompareAny((PetscObject)subsolver, &compatible, TAOSHELL, TAOBNCG, TAOBQNLS, TAOBQNKLS, TAOBQNKTR, TAOBQNKTL, ""));
1263c859ba3SBarry Smith     PetscCheck(compatible,PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "Subsolver must be a bound-constrained first-order method");
127661095bbSAlp Dener   } else {
1289566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompareAny((PetscObject)subsolver, &compatible, TAOSHELL, TAOCG, TAOLMVM, TAOBNCG, TAOBQNLS, TAOBQNKLS, TAOBQNKTR, TAOBQNKTL, ""));
1293c859ba3SBarry Smith     PetscCheck(compatible,PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "Subsolver must be a first-order method");
130661095bbSAlp Dener   }
1313c859ba3SBarry Smith   PetscCheck(compatible,PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "Subsolver must be a first-order method");
1329566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)subsolver));
1339566063dSJacob Faibussowitsch   PetscCall(TaoDestroy(&auglag->subsolver));
134661095bbSAlp Dener   auglag->subsolver = subsolver;
135661095bbSAlp Dener   if (tao->setupcalled) {
1369566063dSJacob Faibussowitsch     PetscCall(TaoSetSolution(auglag->subsolver, auglag->P));
1379566063dSJacob Faibussowitsch     PetscCall(TaoSetObjective(auglag->subsolver, TaoALMMSubsolverObjective_Private, (void*)auglag));
1389566063dSJacob Faibussowitsch     PetscCall(TaoSetObjectiveAndGradient(auglag->subsolver, NULL, TaoALMMSubsolverObjectiveAndGradient_Private, (void*)auglag));
1399566063dSJacob Faibussowitsch     PetscCall(TaoSetVariableBounds(auglag->subsolver, auglag->PL, auglag->PU));
140661095bbSAlp Dener   }
141661095bbSAlp Dener   PetscFunctionReturn(0);
142661095bbSAlp Dener }
143661095bbSAlp Dener 
144661095bbSAlp Dener /*@
145661095bbSAlp Dener    TaoALMMGetMultipliers - Retreive a pointer to the Lagrange multipliers.
146661095bbSAlp Dener 
147661095bbSAlp Dener    Input Parameters:
148661095bbSAlp Dener .  tao - the Tao context for the TAOALMM solver
149661095bbSAlp Dener 
150661095bbSAlp Dener    Output Parameters:
151661095bbSAlp Dener .  Y - vector of Lagrange multipliers
152661095bbSAlp Dener 
153661095bbSAlp Dener    Level: advanced
154661095bbSAlp Dener 
155661095bbSAlp Dener    Notes:
156661095bbSAlp Dener    For problems with both equality and inequality constraints,
157661095bbSAlp Dener    the multipliers are combined together as Y = (Ye, Yi). Users
158661095bbSAlp Dener    can recover copies of the subcomponents using index sets
159661095bbSAlp Dener    provided by TaoALMMGetDualIS() and use VecGetSubVector().
160661095bbSAlp Dener 
161661095bbSAlp Dener .seealso: TAOALMM, TaoALMMSetMultipliers(), TaoALMMGetDualIS()
162661095bbSAlp Dener @*/
163661095bbSAlp Dener PetscErrorCode TaoALMMGetMultipliers(Tao tao, Vec *Y)
164661095bbSAlp Dener {
165661095bbSAlp Dener   PetscFunctionBegin;
166661095bbSAlp Dener   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
167661095bbSAlp Dener   PetscValidPointer(Y, 2);
168*cac4c232SBarry Smith   PetscUseMethod(tao,"TaoALMMGetMultipliers_C",(Tao,Vec *),(tao,Y));
169661095bbSAlp Dener   PetscFunctionReturn(0);
170661095bbSAlp Dener }
171661095bbSAlp Dener 
172661095bbSAlp Dener PetscErrorCode TaoALMMGetMultipliers_Private(Tao tao, Vec *Y)
173661095bbSAlp Dener {
174661095bbSAlp Dener   TAO_ALMM *auglag = (TAO_ALMM*)tao->data;
175661095bbSAlp Dener 
176661095bbSAlp Dener   PetscFunctionBegin;
1773c859ba3SBarry Smith   PetscCheck(tao->setupcalled,PetscObjectComm((PetscObject)tao), PETSC_ERR_ORDER, "TaoSetUp() must be called first for scatters to be constructed");
178661095bbSAlp Dener   *Y = auglag->Y;
179661095bbSAlp Dener   PetscFunctionReturn(0);
180661095bbSAlp Dener }
181661095bbSAlp Dener 
182661095bbSAlp Dener /*@
183661095bbSAlp Dener    TaoALMMSetMultipliers - Set user-defined Lagrange multipliers. The vector type and
184661095bbSAlp Dener                              parallel structure of the given vectormust match equality and
185661095bbSAlp Dener                              inequality constraints. The vector must have a local size equal
186661095bbSAlp Dener                              to the sum of the local sizes for the constraint vectors, and a
187661095bbSAlp Dener                              global size equal to the sum of the global sizes of the constraint
188661095bbSAlp Dener                              vectors.
189661095bbSAlp Dener 
190661095bbSAlp Dener    Input Parameters:
191661095bbSAlp Dener +  tao - the Tao context for the TAOALMM solver
192661095bbSAlp Dener -  Y - vector of Lagrange multipliers
193661095bbSAlp Dener 
194661095bbSAlp Dener    Level: advanced
195661095bbSAlp Dener 
196661095bbSAlp Dener    Notes:
197661095bbSAlp Dener    This routine is only useful if the user wants to change the
198661095bbSAlp Dener    parallel distribution of the combined dual vector in problems that
199661095bbSAlp Dener    feature both equality and inequality constraints. For other tasks,
200661095bbSAlp Dener    it is strongly recommended that the user retreive the dual vector
201661095bbSAlp Dener    created by the solver using TaoALMMGetMultipliers().
202661095bbSAlp Dener 
203661095bbSAlp Dener .seealso: TAOALMM, TaoALMMGetMultipliers()
204661095bbSAlp Dener @*/
205661095bbSAlp Dener PetscErrorCode TaoALMMSetMultipliers(Tao tao, Vec Y)
206661095bbSAlp Dener {
207661095bbSAlp Dener   PetscFunctionBegin;
208661095bbSAlp Dener   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
209661095bbSAlp Dener   PetscValidHeaderSpecific(Y, VEC_CLASSID, 2);
210*cac4c232SBarry Smith   PetscTryMethod(tao,"TaoALMMSetMultipliers_C",(Tao,Vec),(tao,Y));
211661095bbSAlp Dener   PetscFunctionReturn(0);
212661095bbSAlp Dener }
213661095bbSAlp Dener 
214661095bbSAlp Dener PetscErrorCode TaoALMMSetMultipliers_Private(Tao tao, Vec Y)
215661095bbSAlp Dener {
216661095bbSAlp Dener   TAO_ALMM       *auglag = (TAO_ALMM*)tao->data;
217661095bbSAlp Dener   VecType        Ytype;
218661095bbSAlp Dener   PetscInt       Nuser, Neq, Nineq, N;
219661095bbSAlp Dener   PetscBool      same = PETSC_FALSE;
220661095bbSAlp Dener 
221661095bbSAlp Dener   PetscFunctionBegin;
222661095bbSAlp Dener   /* no-op if user provides vector from TaoALMMGetMultipliers() */
223661095bbSAlp Dener   if (Y == auglag->Y) PetscFunctionReturn(0);
224661095bbSAlp Dener   /* make sure vector type is same as equality and inequality constraints */
225661095bbSAlp Dener   if (tao->eq_constrained) {
2269566063dSJacob Faibussowitsch     PetscCall(VecGetType(tao->constraints_equality, &Ytype));
227661095bbSAlp Dener   } else {
2289566063dSJacob Faibussowitsch     PetscCall(VecGetType(tao->constraints_inequality, &Ytype));
229661095bbSAlp Dener   }
2309566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)Y, Ytype, &same));
2313c859ba3SBarry Smith   PetscCheck(same,PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "Given vector for multipliers is not the same type as constraint vectors");
232661095bbSAlp Dener   /* make sure global size matches sum of equality and inequality */
233661095bbSAlp Dener   if (tao->eq_constrained) {
2349566063dSJacob Faibussowitsch     PetscCall(VecGetSize(tao->constraints_equality, &Neq));
235661095bbSAlp Dener   } else {
236661095bbSAlp Dener     Neq = 0;
237661095bbSAlp Dener   }
238661095bbSAlp Dener   if (tao->ineq_constrained) {
2399566063dSJacob Faibussowitsch     PetscCall(VecGetSize(tao->constraints_inequality, &Nineq));
240661095bbSAlp Dener   } else {
241661095bbSAlp Dener     Nineq = 0;
242661095bbSAlp Dener   }
243661095bbSAlp Dener   N = Neq + Nineq;
2449566063dSJacob Faibussowitsch   PetscCall(VecGetSize(Y, &Nuser));
2453c859ba3SBarry Smith   PetscCheck(Nuser == N,PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "Given vector has wrong global size");
246661095bbSAlp Dener   /* if there is only one type of constraint, then we need the local size to match too */
247661095bbSAlp Dener   if (Neq == 0) {
2489566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(tao->constraints_inequality, &Nineq));
2499566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(Y, &Nuser));
2503c859ba3SBarry Smith     PetscCheck(Nuser == Nineq,PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "Given vector has wrong local size");
251661095bbSAlp Dener   }
252661095bbSAlp Dener   if (Nineq == 0) {
2539566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(tao->constraints_equality, &Neq));
2549566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(Y, &Nuser));
2553c859ba3SBarry Smith     PetscCheck(Nuser == Neq,PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_INCOMP, "Given vector has wrong local size");
256661095bbSAlp Dener   }
257661095bbSAlp Dener   /* if we got here, the given vector is compatible so we can replace the current one */
2589566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)Y));
2599566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&auglag->Y));
260661095bbSAlp Dener   auglag->Y = Y;
261661095bbSAlp Dener   /* if there are both types of constraints and the solver has already been set up,
262661095bbSAlp Dener      then we need to regenerate VecScatter objects for the new combined dual vector */
263661095bbSAlp Dener   if (tao->setupcalled && tao->eq_constrained && tao->ineq_constrained) {
2649566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&auglag->C));
2659566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(auglag->Y, &auglag->C));
2669566063dSJacob Faibussowitsch     PetscCall(VecScatterDestroy(&auglag->Yscatter[0]));
2679566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(auglag->Y, auglag->Yis[0], auglag->Ye, NULL, &auglag->Yscatter[0]));
2689566063dSJacob Faibussowitsch     PetscCall(VecScatterDestroy(&auglag->Yscatter[1]));
2699566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(auglag->Y, auglag->Yis[1], auglag->Yi, NULL, &auglag->Yscatter[1]));
270661095bbSAlp Dener   }
271661095bbSAlp Dener   PetscFunctionReturn(0);
272661095bbSAlp Dener }
273661095bbSAlp Dener 
274661095bbSAlp Dener /*@
275661095bbSAlp Dener    TaoALMMGetPrimalIS - Retreive a pointer to the index set that identifies optimization
276661095bbSAlp Dener                         and slack variable components of the subsolver's solution vector.
277661095bbSAlp Dener                         Not valid for problems with only equality constraints.
278661095bbSAlp Dener 
279f899ff85SJose E. Roman    Input Parameter:
280661095bbSAlp Dener .  tao - the Tao context for the TAOALMM solver
281661095bbSAlp Dener 
282661095bbSAlp Dener    Output Parameters:
283661095bbSAlp Dener +  opt_is - index set associated with the optimization variables (NULL if not needed)
284661095bbSAlp Dener -  slack_is - index set associated with the slack variables (NULL if not needed)
285661095bbSAlp Dener 
286661095bbSAlp Dener    Level: advanced
287661095bbSAlp Dener 
288661095bbSAlp Dener .seealso: TAOALMM, TaoALMMGetPrimalVector()
289661095bbSAlp Dener @*/
290661095bbSAlp Dener PetscErrorCode TaoALMMGetPrimalIS(Tao tao, IS *opt_is, IS *slack_is)
291661095bbSAlp Dener {
292661095bbSAlp Dener   PetscFunctionBegin;
293661095bbSAlp Dener   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
294*cac4c232SBarry Smith   PetscUseMethod(tao,"TaoALMMGetPrimalIS_C",(Tao,IS *,IS *),(tao,opt_is,slack_is));
295661095bbSAlp Dener   PetscFunctionReturn(0);
296661095bbSAlp Dener }
297661095bbSAlp Dener 
298661095bbSAlp Dener PetscErrorCode TaoALMMGetPrimalIS_Private(Tao tao, IS *opt_is, IS *slack_is)
299661095bbSAlp Dener {
300661095bbSAlp Dener   TAO_ALMM *auglag = (TAO_ALMM*)tao->data;
301661095bbSAlp Dener 
302661095bbSAlp Dener   PetscFunctionBegin;
3033c859ba3SBarry Smith   PetscCheck(tao->ineq_constrained,PetscObjectComm((PetscObject)tao), PETSC_ERR_ARG_WRONGSTATE, "Primal space has index sets only for inequality constrained problems");
3043c859ba3SBarry Smith   PetscCheck(tao->setupcalled,PetscObjectComm((PetscObject)tao), PETSC_ERR_ORDER, "TaoSetUp() must be called first for index sets to be constructed");
3057721a36fSStefano Zampini   if (opt_is) *opt_is = auglag->Pis[0];
3067721a36fSStefano Zampini   if (slack_is) *slack_is = auglag->Pis[1];
307661095bbSAlp Dener   PetscFunctionReturn(0);
308661095bbSAlp Dener }
309661095bbSAlp Dener 
310661095bbSAlp Dener /*@
311661095bbSAlp Dener    TaoALMMGetDualIS - Retreive a pointer to the index set that identifies equality
312661095bbSAlp Dener                       and inequality constraint components of the dual vector returned
313661095bbSAlp Dener                       by TaoALMMGetMultipliers(). Not valid for problems with only one
314661095bbSAlp Dener                       type of constraint.
315661095bbSAlp Dener 
316f899ff85SJose E. Roman    Input Parameter:
317661095bbSAlp Dener .  tao - the Tao context for the TAOALMM solver
318661095bbSAlp Dener 
319661095bbSAlp Dener    Output Parameters:
320661095bbSAlp Dener +  eq_is - index set associated with the equality constraints (NULL if not needed)
321661095bbSAlp Dener -  ineq_is - index set associated with the inequality constraints (NULL if not needed)
322661095bbSAlp Dener 
323661095bbSAlp Dener    Level: advanced
324661095bbSAlp Dener 
325661095bbSAlp Dener .seealso: TAOALMM, TaoALMMGetMultipliers()
326661095bbSAlp Dener @*/
327661095bbSAlp Dener PetscErrorCode TaoALMMGetDualIS(Tao tao, IS *eq_is, IS *ineq_is)
328661095bbSAlp Dener {
329661095bbSAlp Dener    PetscFunctionBegin;
330661095bbSAlp Dener    PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
331*cac4c232SBarry Smith    PetscUseMethod(tao,"TaoALMMGetDualIS_C",(Tao,IS *,IS *),(tao,eq_is,ineq_is));
332661095bbSAlp Dener    PetscFunctionReturn(0);
333661095bbSAlp Dener }
334661095bbSAlp Dener 
335661095bbSAlp Dener PetscErrorCode TaoALMMGetDualIS_Private(Tao tao, IS *eq_is, IS *ineq_is)
336661095bbSAlp Dener {
337661095bbSAlp Dener   TAO_ALMM *auglag = (TAO_ALMM*)tao->data;
338661095bbSAlp Dener 
339661095bbSAlp Dener   PetscFunctionBegin;
340661095bbSAlp Dener   PetscValidHeaderSpecific(tao, TAO_CLASSID, 1);
3413c859ba3SBarry 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");
3423c859ba3SBarry Smith   PetscCheck(tao->setupcalled,PetscObjectComm((PetscObject)tao), PETSC_ERR_ORDER, "TaoSetUp() must be called first for index sets to be constructed");
3437721a36fSStefano Zampini   if (eq_is) *eq_is = auglag->Yis[0];
3447721a36fSStefano Zampini   if (ineq_is) *ineq_is = auglag->Yis[1];
345661095bbSAlp Dener   PetscFunctionReturn(0);
346661095bbSAlp Dener }
347