xref: /petsc/src/tao/constrained/impls/ipm/ipm.c (revision 1a1499c8e13c12f02cf4c59cfd6b0cfcce01ae9b)
1ba92ff59SBarry Smith #include <petsctaolinesearch.h>
2aaa7dc30SBarry Smith #include <../src/tao/constrained/impls/ipm/ipm.h> /*I "ipm.h" I*/
31522df2eSJason Sarich 
4a7e14dcfSSatish Balay /*
5a7e14dcfSSatish Balay    x,d in R^n
6a7e14dcfSSatish Balay    f in R
7a7e14dcfSSatish Balay    nb = mi + nlb+nub
8a7e14dcfSSatish Balay    s in R^nb is slack vector CI(x) / x-XL / -x+XU
9a7e14dcfSSatish Balay    bin in R^mi (tao->constraints_inequality)
10a7e14dcfSSatish Balay    beq in R^me (tao->constraints_equality)
11a7e14dcfSSatish Balay    lamdai in R^nb (ipmP->lamdai)
12a7e14dcfSSatish Balay    lamdae in R^me (ipmP->lamdae)
13a7e14dcfSSatish Balay    Jeq in R^(me x n) (tao->jacobian_equality)
14a7e14dcfSSatish Balay    Jin in R^(mi x n) (tao->jacobian_inequality)
15a7e14dcfSSatish Balay    Ai in  R^(nb x n) (ipmP->Ai)
16a7e14dcfSSatish Balay    H in R^(n x n) (tao->hessian)
17a7e14dcfSSatish Balay    min f=(1/2)*x'*H*x + d'*x
18a7e14dcfSSatish Balay    s.t.  CE(x) == 0
19a7e14dcfSSatish Balay          CI(x) >= 0
20a7e14dcfSSatish Balay          x >= tao->XL
21a7e14dcfSSatish Balay          -x >= -tao->XU
22a7e14dcfSSatish Balay */
23a7e14dcfSSatish Balay 
24441846f8SBarry Smith static PetscErrorCode IPMComputeKKT(Tao tao);
25441846f8SBarry Smith static PetscErrorCode IPMPushInitialPoint(Tao tao);
26441846f8SBarry Smith static PetscErrorCode IPMEvaluate(Tao tao);
27441846f8SBarry Smith static PetscErrorCode IPMUpdateK(Tao tao);
28441846f8SBarry Smith static PetscErrorCode IPMUpdateAi(Tao tao);
29441846f8SBarry Smith static PetscErrorCode IPMGatherRHS(Tao tao,Vec,Vec,Vec,Vec,Vec);
30441846f8SBarry Smith static PetscErrorCode IPMScatterStep(Tao tao,Vec,Vec,Vec,Vec,Vec);
31441846f8SBarry Smith static PetscErrorCode IPMInitializeBounds(Tao tao);
32a7e14dcfSSatish Balay 
33a7e14dcfSSatish Balay #undef __FUNCT__
34a7e14dcfSSatish Balay #define __FUNCT__ "TaoSolve_IPM"
35441846f8SBarry Smith static PetscErrorCode TaoSolve_IPM(Tao tao)
36a7e14dcfSSatish Balay {
37a7e14dcfSSatish Balay   PetscErrorCode     ierr;
38a7e14dcfSSatish Balay   TAO_IPM            *ipmP = (TAO_IPM*)tao->data;
39e4cb33bbSBarry Smith   TaoConvergedReason reason = TAO_CONTINUE_ITERATING;
40a7e14dcfSSatish Balay   PetscInt           iter = 0,its,i;
41a7e14dcfSSatish Balay   PetscScalar        stepsize=1.0;
42a7e14dcfSSatish Balay   PetscScalar        step_s,step_l,alpha,tau,sigma,phi_target;
43a7e14dcfSSatish Balay 
4447a47007SBarry Smith   PetscFunctionBegin;
45a7e14dcfSSatish Balay   /* Push initial point away from bounds */
46a7e14dcfSSatish Balay   ierr = IPMInitializeBounds(tao);CHKERRQ(ierr);
47a7e14dcfSSatish Balay   ierr = IPMPushInitialPoint(tao);CHKERRQ(ierr);
48a7e14dcfSSatish Balay   ierr = VecCopy(tao->solution,ipmP->rhs_x);CHKERRQ(ierr);
49a7e14dcfSSatish Balay   ierr = IPMEvaluate(tao);CHKERRQ(ierr);
50a7e14dcfSSatish Balay   ierr = IPMComputeKKT(tao);CHKERRQ(ierr);
51a7e14dcfSSatish Balay   ierr = TaoMonitor(tao,iter++,ipmP->kkt_f,ipmP->phi,0.0,1.0,&reason);
52a7e14dcfSSatish Balay 
53a7e14dcfSSatish Balay   while (reason == TAO_CONTINUE_ITERATING) {
54b0026674SJason Sarich     tao->ksp_its=0;
55a7e14dcfSSatish Balay     ierr = IPMUpdateK(tao);CHKERRQ(ierr);
56a7e14dcfSSatish Balay     /*
57a7e14dcfSSatish Balay        rhs.x    = -rd
58a7e14dcfSSatish Balay        rhs.lame = -rpe
59a7e14dcfSSatish Balay        rhs.lami = -rpi
60a7e14dcfSSatish Balay        rhs.com  = -com
61a7e14dcfSSatish Balay     */
62a7e14dcfSSatish Balay 
63a7e14dcfSSatish Balay     ierr = VecCopy(ipmP->rd,ipmP->rhs_x);CHKERRQ(ierr);
64a7e14dcfSSatish Balay     if (ipmP->me > 0) {
65a7e14dcfSSatish Balay       ierr = VecCopy(ipmP->rpe,ipmP->rhs_lamdae);CHKERRQ(ierr);
66a7e14dcfSSatish Balay     }
67a7e14dcfSSatish Balay     if (ipmP->nb > 0) {
68a7e14dcfSSatish Balay       ierr = VecCopy(ipmP->rpi,ipmP->rhs_lamdai);CHKERRQ(ierr);
69a7e14dcfSSatish Balay       ierr = VecCopy(ipmP->complementarity,ipmP->rhs_s);CHKERRQ(ierr);
70a7e14dcfSSatish Balay     }
71050fc7a3SBarry Smith     ierr = IPMGatherRHS(tao,ipmP->bigrhs,ipmP->rhs_x,ipmP->rhs_lamdae,ipmP->rhs_lamdai,ipmP->rhs_s);CHKERRQ(ierr);
72a7e14dcfSSatish Balay     ierr = VecScale(ipmP->bigrhs,-1.0);CHKERRQ(ierr);
73a7e14dcfSSatish Balay 
74a7e14dcfSSatish Balay     /* solve K * step = rhs */
7523ee1639SBarry Smith     ierr = KSPSetOperators(tao->ksp,ipmP->K,ipmP->K);CHKERRQ(ierr);
76a7e14dcfSSatish Balay     ierr = KSPSolve(tao->ksp,ipmP->bigrhs,ipmP->bigstep);CHKERRQ(ierr);
77a7e14dcfSSatish Balay 
78050fc7a3SBarry Smith     ierr = IPMScatterStep(tao,ipmP->bigstep,tao->stepdirection,ipmP->ds,ipmP->dlamdae,ipmP->dlamdai);CHKERRQ(ierr);
79a7e14dcfSSatish Balay     ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr);
80a7e14dcfSSatish Balay     tao->ksp_its += its;
81b0026674SJason Sarich     tao->ksp_tot_its+=its;
82a7e14dcfSSatish Balay      /* Find distance along step direction to closest bound */
83a7e14dcfSSatish Balay     if (ipmP->nb > 0) {
84e270355aSBarry Smith       ierr = VecStepBoundInfo(ipmP->s,ipmP->ds,ipmP->Zero_nb,ipmP->Inf_nb,&step_s,NULL,NULL);CHKERRQ(ierr);
85e270355aSBarry Smith       ierr = VecStepBoundInfo(ipmP->lamdai,ipmP->dlamdai,ipmP->Zero_nb,ipmP->Inf_nb,&step_l,NULL,NULL);CHKERRQ(ierr);
86a7e14dcfSSatish Balay       alpha = PetscMin(step_s,step_l);
87a7e14dcfSSatish Balay       alpha = PetscMin(alpha,1.0);
88a7e14dcfSSatish Balay       ipmP->alpha1 = alpha;
89a7e14dcfSSatish Balay     } else {
90a7e14dcfSSatish Balay       ipmP->alpha1 = alpha = 1.0;
91a7e14dcfSSatish Balay     }
92a7e14dcfSSatish Balay 
93a7e14dcfSSatish Balay     /* x_aff = x + alpha*d */
94a7e14dcfSSatish Balay     ierr = VecCopy(tao->solution,ipmP->save_x);CHKERRQ(ierr);
95a7e14dcfSSatish Balay     if (ipmP->me > 0) {
96a7e14dcfSSatish Balay       ierr = VecCopy(ipmP->lamdae,ipmP->save_lamdae);CHKERRQ(ierr);
97a7e14dcfSSatish Balay     }
98a7e14dcfSSatish Balay     if (ipmP->nb > 0) {
99a7e14dcfSSatish Balay       ierr = VecCopy(ipmP->lamdai,ipmP->save_lamdai);CHKERRQ(ierr);
100a7e14dcfSSatish Balay       ierr = VecCopy(ipmP->s,ipmP->save_s);CHKERRQ(ierr);
101a7e14dcfSSatish Balay     }
102a7e14dcfSSatish Balay 
103a7e14dcfSSatish Balay     ierr = VecAXPY(tao->solution,alpha,tao->stepdirection);CHKERRQ(ierr);
104a7e14dcfSSatish Balay     if (ipmP->me > 0) {
105a7e14dcfSSatish Balay       ierr = VecAXPY(ipmP->lamdae,alpha,ipmP->dlamdae);CHKERRQ(ierr);
106a7e14dcfSSatish Balay     }
107a7e14dcfSSatish Balay     if (ipmP->nb > 0) {
108a7e14dcfSSatish Balay       ierr = VecAXPY(ipmP->lamdai,alpha,ipmP->dlamdai);CHKERRQ(ierr);
109a7e14dcfSSatish Balay       ierr = VecAXPY(ipmP->s,alpha,ipmP->ds);CHKERRQ(ierr);
110a7e14dcfSSatish Balay     }
111a7e14dcfSSatish Balay 
112a7e14dcfSSatish Balay     /* Recompute kkt to find centering parameter sigma = (new_mu/old_mu)^3 */
113a7e14dcfSSatish Balay     if (ipmP->mu == 0.0) {
114a7e14dcfSSatish Balay       sigma = 0.0;
115a7e14dcfSSatish Balay     } else {
116a7e14dcfSSatish Balay       sigma = 1.0/ipmP->mu;
117a7e14dcfSSatish Balay     }
118a7e14dcfSSatish Balay     ierr = IPMComputeKKT(tao);CHKERRQ(ierr);
119a7e14dcfSSatish Balay     sigma *= ipmP->mu;
120a7e14dcfSSatish Balay     sigma*=sigma*sigma;
121a7e14dcfSSatish Balay 
122a7e14dcfSSatish Balay     /* revert kkt info */
123a7e14dcfSSatish Balay     ierr = VecCopy(ipmP->save_x,tao->solution);CHKERRQ(ierr);
124a7e14dcfSSatish Balay     if (ipmP->me > 0) {
125a7e14dcfSSatish Balay       ierr = VecCopy(ipmP->save_lamdae,ipmP->lamdae);CHKERRQ(ierr);
126a7e14dcfSSatish Balay     }
127a7e14dcfSSatish Balay     if (ipmP->nb > 0) {
128a7e14dcfSSatish Balay       ierr = VecCopy(ipmP->save_lamdai,ipmP->lamdai);CHKERRQ(ierr);
129a7e14dcfSSatish Balay       ierr = VecCopy(ipmP->save_s,ipmP->s);CHKERRQ(ierr);
130a7e14dcfSSatish Balay     }
131a7e14dcfSSatish Balay     ierr = IPMComputeKKT(tao);CHKERRQ(ierr);
132a7e14dcfSSatish Balay 
133a7e14dcfSSatish Balay     /* update rhs with new complementarity vector */
134a7e14dcfSSatish Balay     if (ipmP->nb > 0) {
135a7e14dcfSSatish Balay       ierr = VecCopy(ipmP->complementarity,ipmP->rhs_s);CHKERRQ(ierr);
136a7e14dcfSSatish Balay       ierr = VecScale(ipmP->rhs_s,-1.0);CHKERRQ(ierr);
137a7e14dcfSSatish Balay       ierr = VecShift(ipmP->rhs_s,sigma*ipmP->mu);CHKERRQ(ierr);
138a7e14dcfSSatish Balay     }
1396c23d075SBarry Smith     ierr = IPMGatherRHS(tao,ipmP->bigrhs,NULL,NULL,NULL,ipmP->rhs_s);CHKERRQ(ierr);
140a7e14dcfSSatish Balay 
141a7e14dcfSSatish Balay     /* solve K * step = rhs */
14223ee1639SBarry Smith     ierr = KSPSetOperators(tao->ksp,ipmP->K,ipmP->K);CHKERRQ(ierr);
143a7e14dcfSSatish Balay     ierr = KSPSolve(tao->ksp,ipmP->bigrhs,ipmP->bigstep);CHKERRQ(ierr);
144a7e14dcfSSatish Balay 
14547a47007SBarry Smith     ierr = IPMScatterStep(tao,ipmP->bigstep,tao->stepdirection,ipmP->ds,ipmP->dlamdae,ipmP->dlamdai);CHKERRQ(ierr);
146a7e14dcfSSatish Balay     ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr);
147a7e14dcfSSatish Balay     tao->ksp_its += its;
148b0026674SJason Sarich     tao->ksp_tot_its+=its;
149a7e14dcfSSatish Balay     if (ipmP->nb > 0) {
150a7e14dcfSSatish Balay       /* Get max step size and apply frac-to-boundary */
151a7e14dcfSSatish Balay       tau = PetscMax(ipmP->taumin,1.0-ipmP->mu);
152a7e14dcfSSatish Balay       tau = PetscMin(tau,1.0);
153a7e14dcfSSatish Balay       if (tau != 1.0) {
154a7e14dcfSSatish Balay         ierr = VecScale(ipmP->s,tau);CHKERRQ(ierr);
155a7e14dcfSSatish Balay         ierr = VecScale(ipmP->lamdai,tau);CHKERRQ(ierr);
156a7e14dcfSSatish Balay       }
157e270355aSBarry Smith       ierr = VecStepBoundInfo(ipmP->s,ipmP->ds,ipmP->Zero_nb,ipmP->Inf_nb,&step_s,NULL,NULL);CHKERRQ(ierr);
158e270355aSBarry Smith       ierr = VecStepBoundInfo(ipmP->lamdai,ipmP->dlamdai,ipmP->Zero_nb,ipmP->Inf_nb,&step_l,NULL,NULL);CHKERRQ(ierr);
159a7e14dcfSSatish Balay       if (tau != 1.0) {
160a7e14dcfSSatish Balay         ierr = VecCopy(ipmP->save_s,ipmP->s);CHKERRQ(ierr);
161a7e14dcfSSatish Balay         ierr = VecCopy(ipmP->save_lamdai,ipmP->lamdai);CHKERRQ(ierr);
162a7e14dcfSSatish Balay       }
163a7e14dcfSSatish Balay       alpha = PetscMin(step_s,step_l);
164a7e14dcfSSatish Balay       alpha = PetscMin(alpha,1.0);
165a7e14dcfSSatish Balay     } else {
166a7e14dcfSSatish Balay       alpha = 1.0;
167a7e14dcfSSatish Balay     }
168a7e14dcfSSatish Balay     ipmP->alpha2 = alpha;
169a7e14dcfSSatish Balay     /* TODO make phi_target meaningful */
170a7e14dcfSSatish Balay     phi_target = ipmP->dec * ipmP->phi;
171a7e14dcfSSatish Balay     for (i=0; i<11;i++) {
172a7e14dcfSSatish Balay       ierr = VecAXPY(tao->solution,alpha,tao->stepdirection);CHKERRQ(ierr);
173a7e14dcfSSatish Balay       if (ipmP->nb > 0) {
174a7e14dcfSSatish Balay         ierr = VecAXPY(ipmP->s,alpha,ipmP->ds);CHKERRQ(ierr);
175a7e14dcfSSatish Balay         ierr = VecAXPY(ipmP->lamdai,alpha,ipmP->dlamdai);CHKERRQ(ierr);
176a7e14dcfSSatish Balay       }
177a7e14dcfSSatish Balay       if (ipmP->me > 0) {
178a7e14dcfSSatish Balay         ierr = VecAXPY(ipmP->lamdae,alpha,ipmP->dlamdae);CHKERRQ(ierr);
179a7e14dcfSSatish Balay       }
180a7e14dcfSSatish Balay 
181a7e14dcfSSatish Balay       /* update dual variables */
182a7e14dcfSSatish Balay       if (ipmP->me > 0) {
183a7e14dcfSSatish Balay         ierr = VecCopy(ipmP->lamdae,tao->DE);CHKERRQ(ierr);
184a7e14dcfSSatish Balay       }
185a7e14dcfSSatish Balay 
186a7e14dcfSSatish Balay       ierr = IPMEvaluate(tao);CHKERRQ(ierr);
187a7e14dcfSSatish Balay       ierr = IPMComputeKKT(tao);CHKERRQ(ierr);
188a7e14dcfSSatish Balay       if (ipmP->phi <= phi_target) break;
189a7e14dcfSSatish Balay       alpha /= 2.0;
190a7e14dcfSSatish Balay     }
191a7e14dcfSSatish Balay 
19247a47007SBarry Smith     ierr = TaoMonitor(tao,iter,ipmP->kkt_f,ipmP->phi,0.0,stepsize,&reason);CHKERRQ(ierr);
193a7e14dcfSSatish Balay     iter++;
194a7e14dcfSSatish Balay   }
195a7e14dcfSSatish Balay   PetscFunctionReturn(0);
196a7e14dcfSSatish Balay }
197a7e14dcfSSatish Balay 
198a7e14dcfSSatish Balay #undef __FUNCT__
199a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetup_IPM"
200441846f8SBarry Smith static PetscErrorCode TaoSetup_IPM(Tao tao)
201a7e14dcfSSatish Balay {
202a7e14dcfSSatish Balay   TAO_IPM        *ipmP = (TAO_IPM*)tao->data;
203a7e14dcfSSatish Balay   PetscErrorCode ierr;
204a7e14dcfSSatish Balay 
205a7e14dcfSSatish Balay   PetscFunctionBegin;
206a7e14dcfSSatish Balay   ipmP->nb = ipmP->mi = ipmP->me = 0;
207a7e14dcfSSatish Balay   ipmP->K=0;
208a7e14dcfSSatish Balay   ierr = VecGetSize(tao->solution,&ipmP->n);CHKERRQ(ierr);
209a7e14dcfSSatish Balay   if (!tao->gradient) {
210a7e14dcfSSatish Balay     ierr = VecDuplicate(tao->solution, &tao->gradient);CHKERRQ(ierr);
211a7e14dcfSSatish Balay     ierr = VecDuplicate(tao->solution, &tao->stepdirection);CHKERRQ(ierr);
212a7e14dcfSSatish Balay     ierr = VecDuplicate(tao->solution, &ipmP->rd);CHKERRQ(ierr);
213a7e14dcfSSatish Balay     ierr = VecDuplicate(tao->solution, &ipmP->rhs_x);CHKERRQ(ierr);
214a7e14dcfSSatish Balay     ierr = VecDuplicate(tao->solution, &ipmP->work);CHKERRQ(ierr);
215a7e14dcfSSatish Balay     ierr = VecDuplicate(tao->solution, &ipmP->save_x);CHKERRQ(ierr);
216a7e14dcfSSatish Balay   }
217a7e14dcfSSatish Balay   if (tao->constraints_equality) {
218a7e14dcfSSatish Balay     ierr = VecGetSize(tao->constraints_equality,&ipmP->me);CHKERRQ(ierr);
219a7e14dcfSSatish Balay     ierr = VecDuplicate(tao->constraints_equality,&ipmP->lamdae);CHKERRQ(ierr);
220a7e14dcfSSatish Balay     ierr = VecDuplicate(tao->constraints_equality,&ipmP->dlamdae);CHKERRQ(ierr);
221a7e14dcfSSatish Balay     ierr = VecDuplicate(tao->constraints_equality,&ipmP->rhs_lamdae);CHKERRQ(ierr);
222a7e14dcfSSatish Balay     ierr = VecDuplicate(tao->constraints_equality,&ipmP->save_lamdae);CHKERRQ(ierr);
223a7e14dcfSSatish Balay     ierr = VecDuplicate(tao->constraints_equality,&ipmP->rpe);CHKERRQ(ierr);
224a7e14dcfSSatish Balay     ierr = VecDuplicate(tao->constraints_equality,&tao->DE);CHKERRQ(ierr);
225a7e14dcfSSatish Balay   }
226a7e14dcfSSatish Balay   if (tao->constraints_inequality) {
227a7e14dcfSSatish Balay     ierr = VecDuplicate(tao->constraints_inequality,&tao->DI);CHKERRQ(ierr);
228a7e14dcfSSatish Balay   }
229a7e14dcfSSatish Balay   PetscFunctionReturn(0);
230a7e14dcfSSatish Balay }
231a7e14dcfSSatish Balay 
232a7e14dcfSSatish Balay #undef __FUNCT__
233a7e14dcfSSatish Balay #define __FUNCT__ "IPMInitializeBounds"
234441846f8SBarry Smith static PetscErrorCode IPMInitializeBounds(Tao tao)
235a7e14dcfSSatish Balay {
236a7e14dcfSSatish Balay   TAO_IPM        *ipmP = (TAO_IPM*)tao->data;
237a7e14dcfSSatish Balay   Vec            xtmp;
238a7e14dcfSSatish Balay   PetscInt       xstart,xend;
239a7e14dcfSSatish Balay   PetscInt       ucstart,ucend; /* user ci */
240a7e14dcfSSatish Balay   PetscInt       ucestart,uceend; /* user ce */
241a7e14dcfSSatish Balay   PetscInt       sstart,send;
242a7e14dcfSSatish Balay   PetscInt       bigsize;
243a7e14dcfSSatish Balay   PetscInt       i,counter,nloc;
244a7e14dcfSSatish Balay   PetscInt       *cind,*xind,*ucind,*uceind,*stepind;
245a7e14dcfSSatish Balay   VecType        vtype;
246a7e14dcfSSatish Balay   const PetscInt *xli,*xui;
247a7e14dcfSSatish Balay   PetscInt       xl_offset,xu_offset;
248a7e14dcfSSatish Balay   IS             bigxl,bigxu,isuc,isc,isx,sis,is1;
249a7e14dcfSSatish Balay   PetscErrorCode ierr;
250a7e14dcfSSatish Balay   MPI_Comm       comm;
25147a47007SBarry Smith 
252a7e14dcfSSatish Balay   PetscFunctionBegin;
253a7e14dcfSSatish Balay   cind=xind=ucind=uceind=stepind=0;
254a7e14dcfSSatish Balay   ipmP->mi=0;
255a7e14dcfSSatish Balay   ipmP->nxlb=0;
256a7e14dcfSSatish Balay   ipmP->nxub=0;
257a7e14dcfSSatish Balay   ipmP->nb=0;
258a7e14dcfSSatish Balay   ipmP->nslack=0;
259a7e14dcfSSatish Balay 
260a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&xtmp);CHKERRQ(ierr);
261a7e14dcfSSatish Balay   if (!tao->XL && !tao->XU && tao->ops->computebounds) {
262a7e14dcfSSatish Balay     ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr);
263a7e14dcfSSatish Balay   }
264a7e14dcfSSatish Balay   if (tao->XL) {
265e270355aSBarry Smith     ierr = VecSet(xtmp,PETSC_NINFINITY);CHKERRQ(ierr);
266a7e14dcfSSatish Balay     ierr = VecWhichGreaterThan(tao->XL,xtmp,&ipmP->isxl);CHKERRQ(ierr);
267a7e14dcfSSatish Balay     ierr = ISGetSize(ipmP->isxl,&ipmP->nxlb);CHKERRQ(ierr);
268a7e14dcfSSatish Balay   } else {
269a7e14dcfSSatish Balay     ipmP->nxlb=0;
270a7e14dcfSSatish Balay   }
271a7e14dcfSSatish Balay   if (tao->XU) {
272e270355aSBarry Smith     ierr = VecSet(xtmp,PETSC_INFINITY);CHKERRQ(ierr);
273a7e14dcfSSatish Balay     ierr = VecWhichLessThan(tao->XU,xtmp,&ipmP->isxu);CHKERRQ(ierr);
274a7e14dcfSSatish Balay     ierr = ISGetSize(ipmP->isxu,&ipmP->nxub);CHKERRQ(ierr);
275a7e14dcfSSatish Balay   } else {
276a7e14dcfSSatish Balay     ipmP->nxub=0;
277a7e14dcfSSatish Balay   }
278a7e14dcfSSatish Balay   ierr = VecDestroy(&xtmp);CHKERRQ(ierr);
279a7e14dcfSSatish Balay   if (tao->constraints_inequality) {
280a7e14dcfSSatish Balay     ierr = VecGetSize(tao->constraints_inequality,&ipmP->mi);CHKERRQ(ierr);
281a7e14dcfSSatish Balay   } else {
282a7e14dcfSSatish Balay     ipmP->mi = 0;
283a7e14dcfSSatish Balay   }
284a7e14dcfSSatish Balay   ipmP->nb = ipmP->nxlb + ipmP->nxub + ipmP->mi;
285a7e14dcfSSatish Balay 
286a7e14dcfSSatish Balay   comm = ((PetscObject)(tao->solution))->comm;
287a7e14dcfSSatish Balay 
288a7e14dcfSSatish Balay   bigsize = ipmP->n+2*ipmP->nb+ipmP->me;
2890e660641SBarry Smith   ierr = PetscMalloc1(bigsize,&stepind);CHKERRQ(ierr);
2900e660641SBarry Smith   ierr = PetscMalloc1(ipmP->n,&xind);CHKERRQ(ierr);
2910e660641SBarry Smith   ierr = PetscMalloc1(ipmP->me,&uceind);CHKERRQ(ierr);
292a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(tao->solution,&xstart,&xend);CHKERRQ(ierr);
293a7e14dcfSSatish Balay 
294a7e14dcfSSatish Balay   if (ipmP->nb > 0) {
295a7e14dcfSSatish Balay     ierr = VecCreate(comm,&ipmP->s);CHKERRQ(ierr);
296a7e14dcfSSatish Balay     ierr = VecSetSizes(ipmP->s,PETSC_DECIDE,ipmP->nb);CHKERRQ(ierr);
297a7e14dcfSSatish Balay     ierr = VecSetFromOptions(ipmP->s);CHKERRQ(ierr);
298a7e14dcfSSatish Balay     ierr = VecDuplicate(ipmP->s,&ipmP->ds);CHKERRQ(ierr);
299a7e14dcfSSatish Balay     ierr = VecDuplicate(ipmP->s,&ipmP->rhs_s);CHKERRQ(ierr);
300a7e14dcfSSatish Balay     ierr = VecDuplicate(ipmP->s,&ipmP->complementarity);CHKERRQ(ierr);
301a7e14dcfSSatish Balay     ierr = VecDuplicate(ipmP->s,&ipmP->ci);CHKERRQ(ierr);
302a7e14dcfSSatish Balay 
303a7e14dcfSSatish Balay     ierr = VecDuplicate(ipmP->s,&ipmP->lamdai);CHKERRQ(ierr);
304a7e14dcfSSatish Balay     ierr = VecDuplicate(ipmP->s,&ipmP->dlamdai);CHKERRQ(ierr);
305a7e14dcfSSatish Balay     ierr = VecDuplicate(ipmP->s,&ipmP->rhs_lamdai);CHKERRQ(ierr);
306a7e14dcfSSatish Balay     ierr = VecDuplicate(ipmP->s,&ipmP->save_lamdai);CHKERRQ(ierr);
307a7e14dcfSSatish Balay 
308a7e14dcfSSatish Balay     ierr = VecDuplicate(ipmP->s,&ipmP->save_s);CHKERRQ(ierr);
309a7e14dcfSSatish Balay     ierr = VecDuplicate(ipmP->s,&ipmP->rpi);CHKERRQ(ierr);
310a7e14dcfSSatish Balay     ierr = VecDuplicate(ipmP->s,&ipmP->Zero_nb);CHKERRQ(ierr);
311a7e14dcfSSatish Balay     ierr = VecSet(ipmP->Zero_nb,0.0);CHKERRQ(ierr);
312a7e14dcfSSatish Balay     ierr = VecDuplicate(ipmP->s,&ipmP->One_nb);CHKERRQ(ierr);
313a7e14dcfSSatish Balay     ierr = VecSet(ipmP->One_nb,1.0);CHKERRQ(ierr);
314a7e14dcfSSatish Balay     ierr = VecDuplicate(ipmP->s,&ipmP->Inf_nb);CHKERRQ(ierr);
315e270355aSBarry Smith     ierr = VecSet(ipmP->Inf_nb,PETSC_INFINITY);CHKERRQ(ierr);
316a7e14dcfSSatish Balay 
3170e660641SBarry Smith     ierr = PetscMalloc1(ipmP->nb,&cind);CHKERRQ(ierr);
3180e660641SBarry Smith     ierr = PetscMalloc1(ipmP->mi,&ucind);CHKERRQ(ierr);
319a7e14dcfSSatish Balay     ierr = VecGetOwnershipRange(ipmP->s,&sstart,&send);CHKERRQ(ierr);
320a7e14dcfSSatish Balay 
321a7e14dcfSSatish Balay     if (ipmP->mi > 0) {
322a7e14dcfSSatish Balay       ierr = VecGetOwnershipRange(tao->constraints_inequality,&ucstart,&ucend);CHKERRQ(ierr);
323a7e14dcfSSatish Balay       counter=0;
324a7e14dcfSSatish Balay       for (i=ucstart;i<ucend;i++) {
325a7e14dcfSSatish Balay         cind[counter++] = i;
326a7e14dcfSSatish Balay       }
327a7e14dcfSSatish Balay       ierr = ISCreateGeneral(comm,counter,cind,PETSC_COPY_VALUES,&isuc);CHKERRQ(ierr);
328a7e14dcfSSatish Balay       ierr = ISCreateGeneral(comm,counter,cind,PETSC_COPY_VALUES,&isc);CHKERRQ(ierr);
329a7e14dcfSSatish Balay       ierr = VecScatterCreate(tao->constraints_inequality,isuc,ipmP->ci,isc,&ipmP->ci_scat);CHKERRQ(ierr);
330a7e14dcfSSatish Balay 
331a7e14dcfSSatish Balay       ierr = ISDestroy(&isuc);
332a7e14dcfSSatish Balay       ierr = ISDestroy(&isc);
333a7e14dcfSSatish Balay     }
334a7e14dcfSSatish Balay     /* need to know how may xbound indices are on each process */
335a7e14dcfSSatish Balay     /* TODO better way */
336a7e14dcfSSatish Balay     if (ipmP->nxlb) {
337a7e14dcfSSatish Balay       ierr = ISAllGather(ipmP->isxl,&bigxl);CHKERRQ(ierr);
338a7e14dcfSSatish Balay       ierr = ISGetIndices(bigxl,&xli);CHKERRQ(ierr);
339a7e14dcfSSatish Balay       /* find offsets for this processor */
340a7e14dcfSSatish Balay       xl_offset = ipmP->mi;
341a7e14dcfSSatish Balay       for (i=0;i<ipmP->nxlb;i++) {
342a7e14dcfSSatish Balay         if (xli[i] < xstart) {
343a7e14dcfSSatish Balay           xl_offset++;
344a7e14dcfSSatish Balay         } else break;
345a7e14dcfSSatish Balay       }
346a7e14dcfSSatish Balay       ierr = ISRestoreIndices(bigxl,&xli);CHKERRQ(ierr);
347a7e14dcfSSatish Balay 
348a7e14dcfSSatish Balay       ierr = ISGetIndices(ipmP->isxl,&xli);CHKERRQ(ierr);
349a7e14dcfSSatish Balay       ierr = ISGetLocalSize(ipmP->isxl,&nloc);CHKERRQ(ierr);
350a7e14dcfSSatish Balay       for (i=0;i<nloc;i++) {
351a7e14dcfSSatish Balay         xind[i] = xli[i];
352a7e14dcfSSatish Balay         cind[i] = xl_offset+i;
353a7e14dcfSSatish Balay       }
354a7e14dcfSSatish Balay 
355a7e14dcfSSatish Balay       ierr = ISCreateGeneral(comm,nloc,xind,PETSC_COPY_VALUES,&isx);CHKERRQ(ierr);
356a7e14dcfSSatish Balay       ierr = ISCreateGeneral(comm,nloc,cind,PETSC_COPY_VALUES,&isc);CHKERRQ(ierr);
357a7e14dcfSSatish Balay       ierr = VecScatterCreate(tao->XL,isx,ipmP->ci,isc,&ipmP->xl_scat);CHKERRQ(ierr);
358a7e14dcfSSatish Balay       ierr = ISDestroy(&isx);CHKERRQ(ierr);
359a7e14dcfSSatish Balay       ierr = ISDestroy(&isc);CHKERRQ(ierr);
360a7e14dcfSSatish Balay       ierr = ISDestroy(&bigxl);CHKERRQ(ierr);
361a7e14dcfSSatish Balay     }
362a7e14dcfSSatish Balay 
363a7e14dcfSSatish Balay     if (ipmP->nxub) {
364a7e14dcfSSatish Balay       ierr = ISAllGather(ipmP->isxu,&bigxu);CHKERRQ(ierr);
365a7e14dcfSSatish Balay       ierr = ISGetIndices(bigxu,&xui);CHKERRQ(ierr);
366a7e14dcfSSatish Balay       /* find offsets for this processor */
367a7e14dcfSSatish Balay       xu_offset = ipmP->mi + ipmP->nxlb;
368a7e14dcfSSatish Balay       for (i=0;i<ipmP->nxub;i++) {
369a7e14dcfSSatish Balay         if (xui[i] < xstart) {
370a7e14dcfSSatish Balay           xu_offset++;
371a7e14dcfSSatish Balay         } else break;
372a7e14dcfSSatish Balay       }
373a7e14dcfSSatish Balay       ierr = ISRestoreIndices(bigxu,&xui);CHKERRQ(ierr);
374a7e14dcfSSatish Balay 
375a7e14dcfSSatish Balay       ierr = ISGetIndices(ipmP->isxu,&xui);CHKERRQ(ierr);
376a7e14dcfSSatish Balay       ierr = ISGetLocalSize(ipmP->isxu,&nloc);CHKERRQ(ierr);
377a7e14dcfSSatish Balay       for (i=0;i<nloc;i++) {
378a7e14dcfSSatish Balay         xind[i] = xui[i];
379a7e14dcfSSatish Balay         cind[i] = xu_offset+i;
380a7e14dcfSSatish Balay       }
381a7e14dcfSSatish Balay 
382a7e14dcfSSatish Balay       ierr = ISCreateGeneral(comm,nloc,xind,PETSC_COPY_VALUES,&isx);CHKERRQ(ierr);
383a7e14dcfSSatish Balay       ierr = ISCreateGeneral(comm,nloc,cind,PETSC_COPY_VALUES,&isc);CHKERRQ(ierr);
384a7e14dcfSSatish Balay       ierr = VecScatterCreate(tao->XU,isx,ipmP->ci,isc,&ipmP->xu_scat);CHKERRQ(ierr);
385a7e14dcfSSatish Balay       ierr = ISDestroy(&isx);CHKERRQ(ierr);
386a7e14dcfSSatish Balay       ierr = ISDestroy(&isc);CHKERRQ(ierr);
387a7e14dcfSSatish Balay       ierr = ISDestroy(&bigxu);CHKERRQ(ierr);
388a7e14dcfSSatish Balay     }
389a7e14dcfSSatish Balay   }
390a7e14dcfSSatish Balay   ierr = VecCreate(comm,&ipmP->bigrhs);CHKERRQ(ierr);
391a7e14dcfSSatish Balay   ierr = VecGetType(tao->solution,&vtype);CHKERRQ(ierr);
392a7e14dcfSSatish Balay   ierr = VecSetType(ipmP->bigrhs,vtype);CHKERRQ(ierr);
393a7e14dcfSSatish Balay   ierr = VecSetSizes(ipmP->bigrhs,PETSC_DECIDE,bigsize);CHKERRQ(ierr);
394a7e14dcfSSatish Balay   ierr = VecSetFromOptions(ipmP->bigrhs);CHKERRQ(ierr);
395a7e14dcfSSatish Balay   ierr = VecDuplicate(ipmP->bigrhs,&ipmP->bigstep);CHKERRQ(ierr);
396a7e14dcfSSatish Balay 
397a7e14dcfSSatish Balay   /* create scatters for step->x and x->rhs */
398a7e14dcfSSatish Balay   for (i=xstart;i<xend;i++) {
399a7e14dcfSSatish Balay     stepind[i-xstart] = i;
400a7e14dcfSSatish Balay     xind[i-xstart] = i;
401a7e14dcfSSatish Balay   }
402a7e14dcfSSatish Balay   ierr = ISCreateGeneral(comm,xend-xstart,stepind,PETSC_COPY_VALUES,&sis);CHKERRQ(ierr);
403a7e14dcfSSatish Balay   ierr = ISCreateGeneral(comm,xend-xstart,xind,PETSC_COPY_VALUES,&is1);CHKERRQ(ierr);
404a7e14dcfSSatish Balay   ierr = VecScatterCreate(ipmP->bigstep,sis,tao->solution,is1,&ipmP->step1);CHKERRQ(ierr);
405a7e14dcfSSatish Balay   ierr = VecScatterCreate(tao->solution,is1,ipmP->bigrhs,sis,&ipmP->rhs1);CHKERRQ(ierr);
406a7e14dcfSSatish Balay   ierr = ISDestroy(&sis);CHKERRQ(ierr);
407a7e14dcfSSatish Balay   ierr = ISDestroy(&is1);CHKERRQ(ierr);
408a7e14dcfSSatish Balay 
409a7e14dcfSSatish Balay   if (ipmP->nb > 0) {
410a7e14dcfSSatish Balay     for (i=sstart;i<send;i++) {
411a7e14dcfSSatish Balay       stepind[i-sstart] = i+ipmP->n;
412a7e14dcfSSatish Balay       cind[i-sstart] = i;
413a7e14dcfSSatish Balay     }
414a7e14dcfSSatish Balay     ierr = ISCreateGeneral(comm,send-sstart,stepind,PETSC_COPY_VALUES,&sis);CHKERRQ(ierr);
415a7e14dcfSSatish Balay     ierr = ISCreateGeneral(comm,send-sstart,cind,PETSC_COPY_VALUES,&is1);CHKERRQ(ierr);
416a7e14dcfSSatish Balay     ierr = VecScatterCreate(ipmP->bigstep,sis,ipmP->s,is1,&ipmP->step2);CHKERRQ(ierr);
417a7e14dcfSSatish Balay     ierr = ISDestroy(&sis);CHKERRQ(ierr);
418a7e14dcfSSatish Balay 
419a7e14dcfSSatish Balay     for (i=sstart;i<send;i++) {
420a7e14dcfSSatish Balay       stepind[i-sstart] = i+ipmP->n+ipmP->me;
421a7e14dcfSSatish Balay       cind[i-sstart] = i;
422a7e14dcfSSatish Balay     }
423a7e14dcfSSatish Balay     ierr = ISCreateGeneral(comm,send-sstart,stepind,PETSC_COPY_VALUES,&sis);CHKERRQ(ierr);
424a7e14dcfSSatish Balay     ierr = VecScatterCreate(ipmP->s,is1,ipmP->bigrhs,sis,&ipmP->rhs3);CHKERRQ(ierr);
425a7e14dcfSSatish Balay     ierr = ISDestroy(&sis);CHKERRQ(ierr);
426a7e14dcfSSatish Balay     ierr = ISDestroy(&is1);CHKERRQ(ierr);
427a7e14dcfSSatish Balay   }
428a7e14dcfSSatish Balay 
429a7e14dcfSSatish Balay   if (ipmP->me > 0) {
430a7e14dcfSSatish Balay     ierr = VecGetOwnershipRange(tao->constraints_equality,&ucestart,&uceend);CHKERRQ(ierr);
431a7e14dcfSSatish Balay     for (i=ucestart;i<uceend;i++) {
432a7e14dcfSSatish Balay       stepind[i-ucestart] = i + ipmP->n+ipmP->nb;
433a7e14dcfSSatish Balay       uceind[i-ucestart] = i;
434a7e14dcfSSatish Balay     }
435a7e14dcfSSatish Balay 
436a7e14dcfSSatish Balay     ierr = ISCreateGeneral(comm,uceend-ucestart,stepind,PETSC_COPY_VALUES,&sis);CHKERRQ(ierr);
437a7e14dcfSSatish Balay     ierr = ISCreateGeneral(comm,uceend-ucestart,uceind,PETSC_COPY_VALUES,&is1);CHKERRQ(ierr);
438a7e14dcfSSatish Balay     ierr = VecScatterCreate(ipmP->bigstep,sis,tao->constraints_equality,is1,&ipmP->step3);CHKERRQ(ierr);
439a7e14dcfSSatish Balay     ierr = ISDestroy(&sis);CHKERRQ(ierr);
440a7e14dcfSSatish Balay 
441a7e14dcfSSatish Balay     for (i=ucestart;i<uceend;i++) {
442a7e14dcfSSatish Balay       stepind[i-ucestart] = i + ipmP->n;
443a7e14dcfSSatish Balay     }
444a7e14dcfSSatish Balay 
445a7e14dcfSSatish Balay     ierr = ISCreateGeneral(comm,uceend-ucestart,stepind,PETSC_COPY_VALUES,&sis);CHKERRQ(ierr);
446a7e14dcfSSatish Balay     ierr = VecScatterCreate(tao->constraints_equality,is1,ipmP->bigrhs,sis,&ipmP->rhs2);CHKERRQ(ierr);
447a7e14dcfSSatish Balay     ierr = ISDestroy(&sis);CHKERRQ(ierr);
448a7e14dcfSSatish Balay     ierr = ISDestroy(&is1);CHKERRQ(ierr);
449a7e14dcfSSatish Balay   }
450a7e14dcfSSatish Balay 
451a7e14dcfSSatish Balay   if (ipmP->nb > 0) {
452a7e14dcfSSatish Balay     for (i=sstart;i<send;i++) {
453a7e14dcfSSatish Balay       stepind[i-sstart] = i + ipmP->n + ipmP->nb + ipmP->me;
454a7e14dcfSSatish Balay       cind[i-sstart] = i;
455a7e14dcfSSatish Balay     }
456a7e14dcfSSatish Balay     ierr = ISCreateGeneral(comm,send-sstart,cind,PETSC_COPY_VALUES,&is1);
457a7e14dcfSSatish Balay     ierr = ISCreateGeneral(comm,send-sstart,stepind,PETSC_COPY_VALUES,&sis);CHKERRQ(ierr);
458a7e14dcfSSatish Balay     ierr = VecScatterCreate(ipmP->bigstep,sis,ipmP->s,is1,&ipmP->step4);CHKERRQ(ierr);
459a7e14dcfSSatish Balay     ierr = VecScatterCreate(ipmP->s,is1,ipmP->bigrhs,sis,&ipmP->rhs4);CHKERRQ(ierr);
460a7e14dcfSSatish Balay     ierr = ISDestroy(&sis);CHKERRQ(ierr);
461a7e14dcfSSatish Balay     ierr = ISDestroy(&is1);CHKERRQ(ierr);
462a7e14dcfSSatish Balay   }
463a7e14dcfSSatish Balay 
464a7e14dcfSSatish Balay   ierr = PetscFree(stepind);CHKERRQ(ierr);
465a7e14dcfSSatish Balay   ierr = PetscFree(cind);CHKERRQ(ierr);
466a7e14dcfSSatish Balay   ierr = PetscFree(ucind);CHKERRQ(ierr);
467a7e14dcfSSatish Balay   ierr = PetscFree(uceind);CHKERRQ(ierr);
468a7e14dcfSSatish Balay   ierr = PetscFree(xind);CHKERRQ(ierr);
469a7e14dcfSSatish Balay   PetscFunctionReturn(0);
470a7e14dcfSSatish Balay }
471a7e14dcfSSatish Balay 
472a7e14dcfSSatish Balay #undef __FUNCT__
473a7e14dcfSSatish Balay #define __FUNCT__ "TaoDestroy_IPM"
474441846f8SBarry Smith static PetscErrorCode TaoDestroy_IPM(Tao tao)
475a7e14dcfSSatish Balay {
476a7e14dcfSSatish Balay   TAO_IPM        *ipmP = (TAO_IPM*)tao->data;
477a7e14dcfSSatish Balay   PetscErrorCode ierr;
47847a47007SBarry Smith 
479a7e14dcfSSatish Balay   PetscFunctionBegin;
480a7e14dcfSSatish Balay   ierr = VecDestroy(&ipmP->rd);CHKERRQ(ierr);
481a7e14dcfSSatish Balay   ierr = VecDestroy(&ipmP->rpe);CHKERRQ(ierr);
482a7e14dcfSSatish Balay   ierr = VecDestroy(&ipmP->rpi);CHKERRQ(ierr);
483a7e14dcfSSatish Balay   ierr = VecDestroy(&ipmP->work);CHKERRQ(ierr);
484a7e14dcfSSatish Balay   ierr = VecDestroy(&ipmP->lamdae);CHKERRQ(ierr);
485a7e14dcfSSatish Balay   ierr = VecDestroy(&ipmP->lamdai);CHKERRQ(ierr);
486a7e14dcfSSatish Balay   ierr = VecDestroy(&ipmP->s);CHKERRQ(ierr);
487a7e14dcfSSatish Balay   ierr = VecDestroy(&ipmP->ds);CHKERRQ(ierr);
488a7e14dcfSSatish Balay   ierr = VecDestroy(&ipmP->ci);CHKERRQ(ierr);
489a7e14dcfSSatish Balay 
490a7e14dcfSSatish Balay   ierr = VecDestroy(&ipmP->rhs_x);CHKERRQ(ierr);
491a7e14dcfSSatish Balay   ierr = VecDestroy(&ipmP->rhs_lamdae);CHKERRQ(ierr);
492a7e14dcfSSatish Balay   ierr = VecDestroy(&ipmP->rhs_lamdai);CHKERRQ(ierr);
493a7e14dcfSSatish Balay   ierr = VecDestroy(&ipmP->rhs_s);CHKERRQ(ierr);
494a7e14dcfSSatish Balay 
495a7e14dcfSSatish Balay   ierr = VecDestroy(&ipmP->save_x);CHKERRQ(ierr);
496a7e14dcfSSatish Balay   ierr = VecDestroy(&ipmP->save_lamdae);CHKERRQ(ierr);
497a7e14dcfSSatish Balay   ierr = VecDestroy(&ipmP->save_lamdai);CHKERRQ(ierr);
498a7e14dcfSSatish Balay   ierr = VecDestroy(&ipmP->save_s);CHKERRQ(ierr);
499a7e14dcfSSatish Balay 
500a7e14dcfSSatish Balay   ierr = VecScatterDestroy(&ipmP->step1);CHKERRQ(ierr);
501a7e14dcfSSatish Balay   ierr = VecScatterDestroy(&ipmP->step2);CHKERRQ(ierr);
502a7e14dcfSSatish Balay   ierr = VecScatterDestroy(&ipmP->step3);CHKERRQ(ierr);
503a7e14dcfSSatish Balay   ierr = VecScatterDestroy(&ipmP->step4);CHKERRQ(ierr);
504a7e14dcfSSatish Balay 
505a7e14dcfSSatish Balay   ierr = VecScatterDestroy(&ipmP->rhs1);CHKERRQ(ierr);
506a7e14dcfSSatish Balay   ierr = VecScatterDestroy(&ipmP->rhs2);CHKERRQ(ierr);
507a7e14dcfSSatish Balay   ierr = VecScatterDestroy(&ipmP->rhs3);CHKERRQ(ierr);
508a7e14dcfSSatish Balay   ierr = VecScatterDestroy(&ipmP->rhs4);CHKERRQ(ierr);
509a7e14dcfSSatish Balay 
510a7e14dcfSSatish Balay   ierr = VecScatterDestroy(&ipmP->ci_scat);CHKERRQ(ierr);
511a7e14dcfSSatish Balay   ierr = VecScatterDestroy(&ipmP->xl_scat);CHKERRQ(ierr);
512a7e14dcfSSatish Balay   ierr = VecScatterDestroy(&ipmP->xu_scat);CHKERRQ(ierr);
513a7e14dcfSSatish Balay 
514a7e14dcfSSatish Balay   ierr = VecDestroy(&ipmP->dlamdai);CHKERRQ(ierr);
515a7e14dcfSSatish Balay   ierr = VecDestroy(&ipmP->dlamdae);CHKERRQ(ierr);
516a7e14dcfSSatish Balay   ierr = VecDestroy(&ipmP->Zero_nb);CHKERRQ(ierr);
517a7e14dcfSSatish Balay   ierr = VecDestroy(&ipmP->One_nb);CHKERRQ(ierr);
518a7e14dcfSSatish Balay   ierr = VecDestroy(&ipmP->Inf_nb);CHKERRQ(ierr);
519a7e14dcfSSatish Balay   ierr = VecDestroy(&ipmP->complementarity);CHKERRQ(ierr);
520a7e14dcfSSatish Balay 
521a7e14dcfSSatish Balay   ierr = VecDestroy(&ipmP->bigrhs);CHKERRQ(ierr);
522a7e14dcfSSatish Balay   ierr = VecDestroy(&ipmP->bigstep);CHKERRQ(ierr);
523a7e14dcfSSatish Balay   ierr = MatDestroy(&ipmP->Ai);CHKERRQ(ierr);
524a7e14dcfSSatish Balay   ierr = MatDestroy(&ipmP->K);CHKERRQ(ierr);
525a7e14dcfSSatish Balay   ierr = ISDestroy(&ipmP->isxu);CHKERRQ(ierr);
526a7e14dcfSSatish Balay   ierr = ISDestroy(&ipmP->isxl);CHKERRQ(ierr);
527a7e14dcfSSatish Balay   ierr = PetscFree(tao->data);CHKERRQ(ierr);
528a7e14dcfSSatish Balay   PetscFunctionReturn(0);
529a7e14dcfSSatish Balay }
530a7e14dcfSSatish Balay 
531a7e14dcfSSatish Balay #undef __FUNCT__
532a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetFromOptions_IPM"
533*1a1499c8SBarry Smith static PetscErrorCode TaoSetFromOptions_IPM(PetscOptionsObjectType *PetscOptionsObject,Tao tao)
534a7e14dcfSSatish Balay {
535a7e14dcfSSatish Balay   TAO_IPM        *ipmP = (TAO_IPM*)tao->data;
536a7e14dcfSSatish Balay   PetscErrorCode ierr;
53747a47007SBarry Smith 
538a7e14dcfSSatish Balay   PetscFunctionBegin;
539*1a1499c8SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"IPM method for constrained optimization");CHKERRQ(ierr);
54094ae4db5SBarry Smith   ierr = PetscOptionsBool("-tao_ipm_monitorkkt","monitor kkt status",NULL,ipmP->monitorkkt,&ipmP->monitorkkt,NULL);CHKERRQ(ierr);
54194ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ipm_pushs","parameter to push initial slack variables away from bounds",NULL,ipmP->pushs,&ipmP->pushs,NULL);CHKERRQ(ierr);
54294ae4db5SBarry Smith   ierr = PetscOptionsReal("-tao_ipm_pushnu","parameter to push initial (inequality) dual variables away from bounds",NULL,ipmP->pushnu,&ipmP->pushnu,NULL);CHKERRQ(ierr);
543a7e14dcfSSatish Balay   ierr = PetscOptionsTail();CHKERRQ(ierr);
544a7e14dcfSSatish Balay   ierr = KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
545a7e14dcfSSatish Balay   PetscFunctionReturn(0);
546a7e14dcfSSatish Balay }
547a7e14dcfSSatish Balay 
548a7e14dcfSSatish Balay #undef __FUNCT__
549a7e14dcfSSatish Balay #define __FUNCT__ "TaoView_IPM"
550441846f8SBarry Smith static PetscErrorCode TaoView_IPM(Tao tao, PetscViewer viewer)
551a7e14dcfSSatish Balay {
552a7e14dcfSSatish Balay   return 0;
553a7e14dcfSSatish Balay }
554a7e14dcfSSatish Balay 
555a7e14dcfSSatish Balay /* IPMObjectiveAndGradient()
556a7e14dcfSSatish Balay    f = d'x + 0.5 * x' * H * x
557a7e14dcfSSatish Balay    rd = H*x + d + Ae'*lame - Ai'*lami
558a7e14dcfSSatish Balay    rpe = Ae*x - be
559a7e14dcfSSatish Balay    rpi = Ai*x - yi - bi
560a7e14dcfSSatish Balay    mu = yi' * lami/mi;
561a7e14dcfSSatish Balay    com = yi.*lami
562a7e14dcfSSatish Balay 
563a7e14dcfSSatish Balay    phi = ||rd|| + ||rpe|| + ||rpi|| + ||com||
564a7e14dcfSSatish Balay */
565a7e14dcfSSatish Balay /*
566a7e14dcfSSatish Balay #undef __FUNCT__
567a7e14dcfSSatish Balay #define __FUNCT__ "IPMObjective"
568a7e14dcfSSatish Balay static PetscErrorCode IPMObjective(TaoLineSearch ls, Vec X, PetscReal *f, void *tptr)
569a7e14dcfSSatish Balay {
570441846f8SBarry Smith   Tao tao = (Tao)tptr;
571a7e14dcfSSatish Balay   TAO_IPM *ipmP = (TAO_IPM*)tao->data;
572a7e14dcfSSatish Balay   PetscErrorCode ierr;
573a7e14dcfSSatish Balay   PetscFunctionBegin;
574a7e14dcfSSatish Balay   ierr = IPMComputeKKT(tao);CHKERRQ(ierr);
575a7e14dcfSSatish Balay   *f = ipmP->phi;
576a7e14dcfSSatish Balay   PetscFunctionReturn(0);
577a7e14dcfSSatish Balay }
578a7e14dcfSSatish Balay */
579a7e14dcfSSatish Balay 
580a7e14dcfSSatish Balay /*
581a7e14dcfSSatish Balay    f = d'x + 0.5 * x' * H * x
582a7e14dcfSSatish Balay    rd = H*x + d + Ae'*lame - Ai'*lami
583a7e14dcfSSatish Balay        Ai =   jac_ineq
584a7e14dcfSSatish Balay                I (w/lb)
585a7e14dcfSSatish Balay               -I (w/ub)
586a7e14dcfSSatish Balay 
587a7e14dcfSSatish Balay    rpe = ce
588a7e14dcfSSatish Balay    rpi = ci - s;
589a7e14dcfSSatish Balay    com = s.*lami
590a7e14dcfSSatish Balay    mu = yi' * lami/mi;
591a7e14dcfSSatish Balay 
592a7e14dcfSSatish Balay    phi = ||rd|| + ||rpe|| + ||rpi|| + ||com||
593a7e14dcfSSatish Balay */
594a7e14dcfSSatish Balay #undef __FUNCT__
595a7e14dcfSSatish Balay #define __FUNCT__ "IPMComputeKKT"
596441846f8SBarry Smith static PetscErrorCode IPMComputeKKT(Tao tao)
597a7e14dcfSSatish Balay {
598a7e14dcfSSatish Balay   TAO_IPM        *ipmP = (TAO_IPM *)tao->data;
599a7e14dcfSSatish Balay   PetscScalar    norm;
600a7e14dcfSSatish Balay   PetscErrorCode ierr;
60147a47007SBarry Smith 
60247a47007SBarry Smith   PetscFunctionBegin;
603a7e14dcfSSatish Balay   ierr = VecCopy(tao->gradient,ipmP->rd);CHKERRQ(ierr);
604a7e14dcfSSatish Balay 
605a7e14dcfSSatish Balay   if (ipmP->me > 0) {
606a7e14dcfSSatish Balay     /* rd = gradient + Ae'*lamdae */
607a7e14dcfSSatish Balay     ierr = MatMultTranspose(tao->jacobian_equality,ipmP->lamdae,ipmP->work);CHKERRQ(ierr);
608a7e14dcfSSatish Balay     ierr = VecAXPY(ipmP->rd, 1.0, ipmP->work);CHKERRQ(ierr);
609a7e14dcfSSatish Balay 
610a7e14dcfSSatish Balay     /* rpe = ce(x) */
611a7e14dcfSSatish Balay     ierr = VecCopy(tao->constraints_equality,ipmP->rpe);CHKERRQ(ierr);
612a7e14dcfSSatish Balay   }
613a7e14dcfSSatish Balay   if (ipmP->nb > 0) {
614a7e14dcfSSatish Balay     /* rd = rd - Ai'*lamdai */
615a7e14dcfSSatish Balay     ierr = MatMultTranspose(ipmP->Ai,ipmP->lamdai,ipmP->work);CHKERRQ(ierr);
616a7e14dcfSSatish Balay     ierr = VecAXPY(ipmP->rd, -1.0, ipmP->work);CHKERRQ(ierr);
6171522df2eSJason Sarich 
618a7e14dcfSSatish Balay     /* rpi = cin - s */
619a7e14dcfSSatish Balay     ierr = VecCopy(ipmP->ci,ipmP->rpi);CHKERRQ(ierr);
620a7e14dcfSSatish Balay     ierr = VecAXPY(ipmP->rpi, -1.0, ipmP->s);CHKERRQ(ierr);
621a7e14dcfSSatish Balay 
622a7e14dcfSSatish Balay     /* com = s .* lami */
623a7e14dcfSSatish Balay     ierr = VecPointwiseMult(ipmP->complementarity, ipmP->s,ipmP->lamdai);CHKERRQ(ierr);
624a7e14dcfSSatish Balay   }
625a7e14dcfSSatish Balay   /* phi = ||rd; rpe; rpi; com|| */
626a7e14dcfSSatish Balay   ierr = VecDot(ipmP->rd,ipmP->rd,&norm);CHKERRQ(ierr);
627a7e14dcfSSatish Balay   ipmP->phi = norm;
628a7e14dcfSSatish Balay   if (ipmP->me > 0 ) {
629a7e14dcfSSatish Balay     ierr = VecDot(ipmP->rpe,ipmP->rpe,&norm);CHKERRQ(ierr);
630a7e14dcfSSatish Balay     ipmP->phi += norm;
631a7e14dcfSSatish Balay   }
632a7e14dcfSSatish Balay   if (ipmP->nb > 0) {
633a7e14dcfSSatish Balay     ierr = VecDot(ipmP->rpi,ipmP->rpi,&norm);CHKERRQ(ierr);
634a7e14dcfSSatish Balay     ipmP->phi += norm;
635a7e14dcfSSatish Balay     ierr = VecDot(ipmP->complementarity,ipmP->complementarity,&norm);CHKERRQ(ierr);
636a7e14dcfSSatish Balay     ipmP->phi += norm;
637a7e14dcfSSatish Balay     /* mu = s'*lami/nb */
638a7e14dcfSSatish Balay     ierr = VecDot(ipmP->s,ipmP->lamdai,&ipmP->mu);CHKERRQ(ierr);
639a7e14dcfSSatish Balay     ipmP->mu /= ipmP->nb;
640a7e14dcfSSatish Balay   } else {
641a7e14dcfSSatish Balay     ipmP->mu = 1.0;
642a7e14dcfSSatish Balay   }
643a7e14dcfSSatish Balay 
644a7e14dcfSSatish Balay   ipmP->phi = PetscSqrtScalar(ipmP->phi);
645a7e14dcfSSatish Balay   PetscFunctionReturn(0);
646a7e14dcfSSatish Balay }
647a7e14dcfSSatish Balay 
648a7e14dcfSSatish Balay #undef __FUNCT__
649a7e14dcfSSatish Balay #define __FUNCT__ "IPMEvaluate"
650a7e14dcfSSatish Balay /* evaluate user info at current point */
651441846f8SBarry Smith PetscErrorCode IPMEvaluate(Tao tao)
652a7e14dcfSSatish Balay {
653a7e14dcfSSatish Balay   TAO_IPM        *ipmP = (TAO_IPM *)tao->data;
654a7e14dcfSSatish Balay   PetscErrorCode ierr;
65547a47007SBarry Smith 
656a7e14dcfSSatish Balay   PetscFunctionBegin;
657a7e14dcfSSatish Balay   ierr = TaoComputeObjectiveAndGradient(tao,tao->solution,&ipmP->kkt_f,tao->gradient);CHKERRQ(ierr);
658ffad9901SBarry Smith   ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
659a7e14dcfSSatish Balay   if (ipmP->me > 0) {
660a7e14dcfSSatish Balay     ierr = TaoComputeEqualityConstraints(tao,tao->solution,tao->constraints_equality);
661ffad9901SBarry Smith     ierr = TaoComputeJacobianEquality(tao,tao->solution,tao->jacobian_equality,tao->jacobian_equality_pre);CHKERRQ(ierr);
662a7e14dcfSSatish Balay   }
663a7e14dcfSSatish Balay   if (ipmP->mi > 0) {
664a7e14dcfSSatish Balay     ierr = TaoComputeInequalityConstraints(tao,tao->solution,tao->constraints_inequality);
665ffad9901SBarry Smith     ierr = TaoComputeJacobianInequality(tao,tao->solution,tao->jacobian_inequality,tao->jacobian_inequality_pre);CHKERRQ(ierr);
666a7e14dcfSSatish Balay   }
667a7e14dcfSSatish Balay   if (ipmP->nb > 0) {
668a7e14dcfSSatish Balay     /* Ai' =   jac_ineq | I (w/lb) | -I (w/ub)  */
669a7e14dcfSSatish Balay     ierr = IPMUpdateAi(tao);CHKERRQ(ierr);
670a7e14dcfSSatish Balay   }
671a7e14dcfSSatish Balay   PetscFunctionReturn(0);
672a7e14dcfSSatish Balay }
673a7e14dcfSSatish Balay 
674a7e14dcfSSatish Balay #undef __FUNCT__
675a7e14dcfSSatish Balay #define __FUNCT__ "IPMPushInitialPoint"
676a7e14dcfSSatish Balay /* Push initial point away from bounds */
677441846f8SBarry Smith PetscErrorCode IPMPushInitialPoint(Tao tao)
678a7e14dcfSSatish Balay {
679a7e14dcfSSatish Balay   TAO_IPM        *ipmP = (TAO_IPM *)tao->data;
680a7e14dcfSSatish Balay   PetscErrorCode ierr;
681a7e14dcfSSatish Balay 
68247a47007SBarry Smith   PetscFunctionBegin;
683a7e14dcfSSatish Balay   ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr);
684a7e14dcfSSatish Balay   if (tao->XL && tao->XU) {
685a7e14dcfSSatish Balay     ierr = VecMedian(tao->XL, tao->solution, tao->XU, tao->solution);CHKERRQ(ierr);
686a7e14dcfSSatish Balay   }
687a7e14dcfSSatish Balay   if (ipmP->nb > 0) {
688a7e14dcfSSatish Balay     ierr = VecSet(ipmP->s,ipmP->pushs);CHKERRQ(ierr);
689a7e14dcfSSatish Balay     ierr = VecSet(ipmP->lamdai,ipmP->pushnu);CHKERRQ(ierr);
690a7e14dcfSSatish Balay     if (ipmP->mi > 0) {
691a7e14dcfSSatish Balay       ierr = VecSet(tao->DI,ipmP->pushnu);CHKERRQ(ierr);
692a7e14dcfSSatish Balay     }
693a7e14dcfSSatish Balay   }
694a7e14dcfSSatish Balay   if (ipmP->me > 0) {
695a7e14dcfSSatish Balay     ierr = VecSet(tao->DE,1.0);CHKERRQ(ierr);
696a7e14dcfSSatish Balay     ierr = VecSet(ipmP->lamdae,1.0);CHKERRQ(ierr);
697a7e14dcfSSatish Balay   }
698a7e14dcfSSatish Balay   PetscFunctionReturn(0);
699a7e14dcfSSatish Balay }
700a7e14dcfSSatish Balay 
701a7e14dcfSSatish Balay #undef __FUNCT__
702a7e14dcfSSatish Balay #define __FUNCT__ "IPMUpdateAi"
703441846f8SBarry Smith PetscErrorCode IPMUpdateAi(Tao tao)
704a7e14dcfSSatish Balay {
705a7e14dcfSSatish Balay   /* Ai =     Ji
706a7e14dcfSSatish Balay               I (w/lb)
707a7e14dcfSSatish Balay              -I (w/ub) */
708a7e14dcfSSatish Balay 
709a7e14dcfSSatish Balay   /* Ci =    user->ci
710a7e14dcfSSatish Balay              Xi - lb (w/lb)
711a7e14dcfSSatish Balay              -Xi + ub (w/ub)  */
712a7e14dcfSSatish Balay 
713a7e14dcfSSatish Balay   TAO_IPM           *ipmP = (TAO_IPM *)tao->data;
714a7e14dcfSSatish Balay   MPI_Comm          comm;
715a7e14dcfSSatish Balay   PetscInt          i;
716a7e14dcfSSatish Balay   PetscScalar       newval;
717a7e14dcfSSatish Balay   PetscInt          newrow,newcol,ncols;
718a7e14dcfSSatish Balay   const PetscScalar *vals;
719a7e14dcfSSatish Balay   const PetscInt    *cols;
720a7e14dcfSSatish Balay   PetscInt          astart,aend,jstart,jend;
721a7e14dcfSSatish Balay   PetscInt          *nonzeros;
722a7e14dcfSSatish Balay   PetscInt          r2,r3,r4;
72347a47007SBarry Smith   PetscMPIInt       size;
724a7e14dcfSSatish Balay   PetscErrorCode    ierr;
725a7e14dcfSSatish Balay 
726a7e14dcfSSatish Balay   PetscFunctionBegin;
727a7e14dcfSSatish Balay   r2 = ipmP->mi;
728a7e14dcfSSatish Balay   r3 = r2 + ipmP->nxlb;
729a7e14dcfSSatish Balay   r4 = r3 + ipmP->nxub;
730a7e14dcfSSatish Balay 
731050fc7a3SBarry Smith   if (!ipmP->nb) PetscFunctionReturn(0);
732a7e14dcfSSatish Balay 
733a7e14dcfSSatish Balay   /* Create Ai matrix if it doesn't exist yet */
734a7e14dcfSSatish Balay   if (!ipmP->Ai) {
735a7e14dcfSSatish Balay     comm = ((PetscObject)(tao->solution))->comm;
7360e660641SBarry Smith     ierr = PetscMalloc1(ipmP->nb,&nonzeros);CHKERRQ(ierr);
73747a47007SBarry Smith     ierr = MPI_Comm_size(comm,&size);
73847a47007SBarry Smith     if (size == 1) {
739a7e14dcfSSatish Balay       for (i=0;i<ipmP->mi;i++) {
7406c23d075SBarry Smith         ierr = MatGetRow(tao->jacobian_inequality,i,&ncols,NULL,NULL);CHKERRQ(ierr);
741a7e14dcfSSatish Balay         nonzeros[i] = ncols;
7426c23d075SBarry Smith         ierr = MatRestoreRow(tao->jacobian_inequality,i,&ncols,NULL,NULL);CHKERRQ(ierr);
743a7e14dcfSSatish Balay       }
744a7e14dcfSSatish Balay       for (i=r2;i<r4;i++) {
745a7e14dcfSSatish Balay         nonzeros[i] = 1;
746a7e14dcfSSatish Balay       }
747a7e14dcfSSatish Balay     }
748a7e14dcfSSatish Balay     ierr = MatCreate(comm,&ipmP->Ai);CHKERRQ(ierr);
749a7e14dcfSSatish Balay     ierr = MatSetType(ipmP->Ai,MATAIJ);CHKERRQ(ierr);
750a7e14dcfSSatish Balay     ierr = MatSetSizes(ipmP->Ai,PETSC_DECIDE,PETSC_DECIDE,ipmP->nb,ipmP->n);CHKERRQ(ierr);
751a7e14dcfSSatish Balay     ierr = MatSetFromOptions(ipmP->Ai);CHKERRQ(ierr);
7526c23d075SBarry Smith     ierr = MatMPIAIJSetPreallocation(ipmP->Ai,ipmP->nb,NULL,ipmP->nb,NULL);
753a7e14dcfSSatish Balay     ierr = MatSeqAIJSetPreallocation(ipmP->Ai,PETSC_DEFAULT,nonzeros);CHKERRQ(ierr);
75447a47007SBarry Smith     if (size ==1) {
755a7e14dcfSSatish Balay       ierr = PetscFree(nonzeros);CHKERRQ(ierr);
756a7e14dcfSSatish Balay     }
757a7e14dcfSSatish Balay   }
758a7e14dcfSSatish Balay 
759a7e14dcfSSatish Balay   /* Copy values from user jacobian to Ai */
760a7e14dcfSSatish Balay   ierr = MatGetOwnershipRange(ipmP->Ai,&astart,&aend);CHKERRQ(ierr);
761a7e14dcfSSatish Balay 
762a7e14dcfSSatish Balay   /* Ai w/lb */
763a7e14dcfSSatish Balay   if (ipmP->mi) {
764a7e14dcfSSatish Balay     ierr = MatZeroEntries(ipmP->Ai);CHKERRQ(ierr);
765a7e14dcfSSatish Balay     ierr = MatGetOwnershipRange(tao->jacobian_inequality,&jstart,&jend);CHKERRQ(ierr);
766a7e14dcfSSatish Balay     for (i=jstart;i<jend;i++) {
767a7e14dcfSSatish Balay       ierr = MatGetRow(tao->jacobian_inequality,i,&ncols,&cols,&vals);CHKERRQ(ierr);
768a7e14dcfSSatish Balay       newrow = i;
769a7e14dcfSSatish Balay       ierr = MatSetValues(ipmP->Ai,1,&newrow,ncols,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
770a7e14dcfSSatish Balay       ierr = MatRestoreRow(tao->jacobian_inequality,i,&ncols,&cols,&vals);CHKERRQ(ierr);
771a7e14dcfSSatish Balay     }
772a7e14dcfSSatish Balay   }
773a7e14dcfSSatish Balay 
774a7e14dcfSSatish Balay   /* I w/ xlb */
775a7e14dcfSSatish Balay   if (ipmP->nxlb) {
776a7e14dcfSSatish Balay     for (i=0;i<ipmP->nxlb;i++) {
777a7e14dcfSSatish Balay       if (i>=astart && i<aend) {
778a7e14dcfSSatish Balay         newrow = i+r2;
779a7e14dcfSSatish Balay         newcol = i;
780a7e14dcfSSatish Balay         newval = 1.0;
781a7e14dcfSSatish Balay         ierr = MatSetValues(ipmP->Ai,1,&newrow,1,&newcol,&newval,INSERT_VALUES);CHKERRQ(ierr);
782a7e14dcfSSatish Balay       }
783a7e14dcfSSatish Balay     }
784a7e14dcfSSatish Balay   }
785a7e14dcfSSatish Balay   if (ipmP->nxub) {
786a7e14dcfSSatish Balay     /* I w/ xub */
787a7e14dcfSSatish Balay     for (i=0;i<ipmP->nxub;i++) {
788a7e14dcfSSatish Balay       if (i>=astart && i<aend) {
789a7e14dcfSSatish Balay       newrow = i+r3;
790a7e14dcfSSatish Balay       newcol = i;
791a7e14dcfSSatish Balay       newval = -1.0;
792a7e14dcfSSatish Balay       ierr = MatSetValues(ipmP->Ai,1,&newrow,1,&newcol,&newval,INSERT_VALUES);CHKERRQ(ierr);
793a7e14dcfSSatish Balay       }
794a7e14dcfSSatish Balay     }
795a7e14dcfSSatish Balay   }
796a7e14dcfSSatish Balay 
797a7e14dcfSSatish Balay   ierr = MatAssemblyBegin(ipmP->Ai,MAT_FINAL_ASSEMBLY);
798a7e14dcfSSatish Balay   ierr = MatAssemblyEnd(ipmP->Ai,MAT_FINAL_ASSEMBLY);
799a7e14dcfSSatish Balay   CHKMEMQ;
800a7e14dcfSSatish Balay 
801a7e14dcfSSatish Balay   ierr = VecSet(ipmP->ci,0.0);CHKERRQ(ierr);
802a7e14dcfSSatish Balay 
803a7e14dcfSSatish Balay   /* user ci */
804a7e14dcfSSatish Balay   if (ipmP->mi > 0) {
805a7e14dcfSSatish Balay     ierr = VecScatterBegin(ipmP->ci_scat,tao->constraints_inequality,ipmP->ci,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
806a7e14dcfSSatish Balay     ierr = VecScatterEnd(ipmP->ci_scat,tao->constraints_inequality,ipmP->ci,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
807a7e14dcfSSatish Balay   }
808a7e14dcfSSatish Balay   if (!ipmP->work){
809a7e14dcfSSatish Balay     VecDuplicate(tao->solution,&ipmP->work);
810a7e14dcfSSatish Balay   }
811a7e14dcfSSatish Balay   ierr = VecCopy(tao->solution,ipmP->work);CHKERRQ(ierr);
812a7e14dcfSSatish Balay   if (tao->XL) {
813a7e14dcfSSatish Balay     ierr = VecAXPY(ipmP->work,-1.0,tao->XL);CHKERRQ(ierr);
814a7e14dcfSSatish Balay 
815a7e14dcfSSatish Balay     /* lower bounds on variables */
816a7e14dcfSSatish Balay     if (ipmP->nxlb > 0) {
817a7e14dcfSSatish Balay       ierr = VecScatterBegin(ipmP->xl_scat,ipmP->work,ipmP->ci,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
818a7e14dcfSSatish Balay       ierr = VecScatterEnd(ipmP->xl_scat,ipmP->work,ipmP->ci,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
819a7e14dcfSSatish Balay     }
820a7e14dcfSSatish Balay   }
821a7e14dcfSSatish Balay   if (tao->XU) {
822a7e14dcfSSatish Balay     /* upper bounds on variables */
823a7e14dcfSSatish Balay     ierr = VecCopy(tao->solution,ipmP->work);CHKERRQ(ierr);
824a7e14dcfSSatish Balay     ierr = VecScale(ipmP->work,-1.0);CHKERRQ(ierr);
825a7e14dcfSSatish Balay     ierr = VecAXPY(ipmP->work,1.0,tao->XU);CHKERRQ(ierr);
826a7e14dcfSSatish Balay     if (ipmP->nxub > 0) {
827a7e14dcfSSatish Balay       ierr = VecScatterBegin(ipmP->xu_scat,ipmP->work,ipmP->ci,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
828a7e14dcfSSatish Balay       ierr = VecScatterEnd(ipmP->xu_scat,ipmP->work,ipmP->ci,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
829a7e14dcfSSatish Balay     }
830a7e14dcfSSatish Balay   }
831a7e14dcfSSatish Balay   PetscFunctionReturn(0);
832a7e14dcfSSatish Balay }
833a7e14dcfSSatish Balay 
834a7e14dcfSSatish Balay #undef __FUNCT__
835a7e14dcfSSatish Balay #define __FUNCT__ "IPMUpdateK"
836a7e14dcfSSatish Balay /* create K = [ Hlag , 0 , Ae', -Ai'];
837a7e14dcfSSatish Balay               [Ae , 0,   0  , 0];
838a7e14dcfSSatish Balay               [Ai ,-I,   0 ,  0];
839a7e14dcfSSatish Balay               [ 0 , S ,  0,   Y ];  */
840441846f8SBarry Smith PetscErrorCode IPMUpdateK(Tao tao)
841a7e14dcfSSatish Balay {
842a7e14dcfSSatish Balay   TAO_IPM         *ipmP = (TAO_IPM *)tao->data;
843a7e14dcfSSatish Balay   MPI_Comm        comm;
84447a47007SBarry Smith   PetscMPIInt     size;
845a7e14dcfSSatish Balay   PetscErrorCode  ierr;
846a7e14dcfSSatish Balay   PetscInt        i,j,row;
847a7e14dcfSSatish Balay   PetscInt        ncols,newcol,newcols[2],newrow;
848a7e14dcfSSatish Balay   const PetscInt  *cols;
849a7e14dcfSSatish Balay   const PetscReal *vals;
8505e081366SBarry Smith   const PetscReal *l,*y;
851a7e14dcfSSatish Balay   PetscReal       *newvals;
852a7e14dcfSSatish Balay   PetscReal       newval;
853a7e14dcfSSatish Balay   PetscInt        subsize;
854a7e14dcfSSatish Balay   const PetscInt  *indices;
855a7e14dcfSSatish Balay   PetscInt        *nonzeros,*d_nonzeros,*o_nonzeros;
856a7e14dcfSSatish Balay   PetscInt        bigsize;
857a7e14dcfSSatish Balay   PetscInt        r1,r2,r3;
858a7e14dcfSSatish Balay   PetscInt        c1,c2,c3;
859a7e14dcfSSatish Balay   PetscInt        klocalsize;
860a7e14dcfSSatish Balay   PetscInt        hstart,hend,kstart,kend;
861a7e14dcfSSatish Balay   PetscInt        aistart,aiend,aestart,aeend;
862a7e14dcfSSatish Balay   PetscInt        sstart,send;
863a7e14dcfSSatish Balay 
86447a47007SBarry Smith   PetscFunctionBegin;
865a7e14dcfSSatish Balay   comm = ((PetscObject)(tao->solution))->comm;
866050fc7a3SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
867a7e14dcfSSatish Balay   ierr = IPMUpdateAi(tao);CHKERRQ(ierr);
8681522df2eSJason Sarich 
869a7e14dcfSSatish Balay   /* allocate workspace */
870a7e14dcfSSatish Balay   subsize = PetscMax(ipmP->n,ipmP->nb);
871a7e14dcfSSatish Balay   subsize = PetscMax(ipmP->me,subsize);
872a7e14dcfSSatish Balay   subsize = PetscMax(2,subsize);
8730e660641SBarry Smith   ierr = PetscMalloc1(subsize,&indices);CHKERRQ(ierr);
8740e660641SBarry Smith   ierr = PetscMalloc1(subsize,&newvals);CHKERRQ(ierr);
875a7e14dcfSSatish Balay 
876a7e14dcfSSatish Balay   r1 = c1 = ipmP->n;
877a7e14dcfSSatish Balay   r2 = r1 + ipmP->me;  c2 = c1 + ipmP->nb;
878a7e14dcfSSatish Balay   r3 = c3 = r2 + ipmP->nb;
879a7e14dcfSSatish Balay 
880a7e14dcfSSatish Balay   bigsize = ipmP->n+2*ipmP->nb+ipmP->me;
881a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(ipmP->bigrhs,&kstart,&kend);CHKERRQ(ierr);
882a7e14dcfSSatish Balay   ierr = MatGetOwnershipRange(tao->hessian,&hstart,&hend);CHKERRQ(ierr);
883a7e14dcfSSatish Balay   klocalsize = kend-kstart;
884a7e14dcfSSatish Balay   if (!ipmP->K) {
88547a47007SBarry Smith     if (size == 1) {
886854ce69bSBarry Smith       ierr = PetscMalloc1(kend-kstart,&nonzeros);CHKERRQ(ierr);
887a7e14dcfSSatish Balay       for (i=0;i<bigsize;i++) {
888a7e14dcfSSatish Balay         if (i<r1) {
8896c23d075SBarry Smith           ierr = MatGetRow(tao->hessian,i,&ncols,NULL,NULL);CHKERRQ(ierr);
890a7e14dcfSSatish Balay           nonzeros[i] = ncols;
8916c23d075SBarry Smith           ierr = MatRestoreRow(tao->hessian,i,&ncols,NULL,NULL);CHKERRQ(ierr);
892a7e14dcfSSatish Balay           nonzeros[i] += ipmP->me+ipmP->nb;
893a7e14dcfSSatish Balay         } else if (i<r2) {
894a7e14dcfSSatish Balay           nonzeros[i-kstart] = ipmP->n;
895a7e14dcfSSatish Balay         } else if (i<r3) {
896a7e14dcfSSatish Balay           nonzeros[i-kstart] = ipmP->n+1;
897a7e14dcfSSatish Balay         } else if (i<bigsize) {
898a7e14dcfSSatish Balay           nonzeros[i-kstart] = 2;
899a7e14dcfSSatish Balay         }
900a7e14dcfSSatish Balay       }
901a7e14dcfSSatish Balay       ierr = MatCreate(comm,&ipmP->K);CHKERRQ(ierr);
902a7e14dcfSSatish Balay       ierr = MatSetType(ipmP->K,MATSEQAIJ);CHKERRQ(ierr);
903a7e14dcfSSatish Balay       ierr = MatSetSizes(ipmP->K,klocalsize,klocalsize,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
904a7e14dcfSSatish Balay       ierr = MatSeqAIJSetPreallocation(ipmP->K,0,nonzeros);CHKERRQ(ierr);
905a7e14dcfSSatish Balay       ierr = MatSetFromOptions(ipmP->K);CHKERRQ(ierr);
906a7e14dcfSSatish Balay       ierr = PetscFree(nonzeros);CHKERRQ(ierr);
907a7e14dcfSSatish Balay     } else {
908854ce69bSBarry Smith       ierr = PetscMalloc1(kend-kstart,&d_nonzeros);CHKERRQ(ierr);
909854ce69bSBarry Smith       ierr = PetscMalloc1(kend-kstart,&o_nonzeros);CHKERRQ(ierr);
910a7e14dcfSSatish Balay       for (i=kstart;i<kend;i++) {
911a7e14dcfSSatish Balay         if (i<r1) {
912a7e14dcfSSatish Balay           /* TODO fix preallocation for mpi mats */
913a7e14dcfSSatish Balay           d_nonzeros[i-kstart] = PetscMin(ipmP->n+ipmP->me+ipmP->nb,kend-kstart);
914a7e14dcfSSatish Balay           o_nonzeros[i-kstart] = PetscMin(ipmP->n+ipmP->me+ipmP->nb,bigsize-(kend-kstart));
915a7e14dcfSSatish Balay         } else if (i<r2) {
916a7e14dcfSSatish Balay           d_nonzeros[i-kstart] = PetscMin(ipmP->n,kend-kstart);
917a7e14dcfSSatish Balay           o_nonzeros[i-kstart] = PetscMin(ipmP->n,bigsize-(kend-kstart));
918a7e14dcfSSatish Balay         } else if (i<r3) {
919a7e14dcfSSatish Balay           d_nonzeros[i-kstart] = PetscMin(ipmP->n+2,kend-kstart);
920a7e14dcfSSatish Balay           o_nonzeros[i-kstart] = PetscMin(ipmP->n+2,bigsize-(kend-kstart));
921a7e14dcfSSatish Balay         } else {
922a7e14dcfSSatish Balay           d_nonzeros[i-kstart] = PetscMin(2,kend-kstart);
923a7e14dcfSSatish Balay           o_nonzeros[i-kstart] = PetscMin(2,bigsize-(kend-kstart));
924a7e14dcfSSatish Balay         }
925a7e14dcfSSatish Balay       }
926a7e14dcfSSatish Balay       ierr = MatCreate(comm,&ipmP->K);CHKERRQ(ierr);
927a7e14dcfSSatish Balay       ierr = MatSetType(ipmP->K,MATMPIAIJ);CHKERRQ(ierr);
928a7e14dcfSSatish Balay       ierr = MatSetSizes(ipmP->K,klocalsize,klocalsize,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
929a7e14dcfSSatish Balay       ierr = MatMPIAIJSetPreallocation(ipmP->K,0,d_nonzeros,0,o_nonzeros);CHKERRQ(ierr);
930a7e14dcfSSatish Balay       ierr = PetscFree(d_nonzeros);CHKERRQ(ierr);
931a7e14dcfSSatish Balay       ierr = PetscFree(o_nonzeros);CHKERRQ(ierr);
932a7e14dcfSSatish Balay       ierr = MatSetFromOptions(ipmP->K);CHKERRQ(ierr);
933a7e14dcfSSatish Balay     }
934a7e14dcfSSatish Balay   }
935a7e14dcfSSatish Balay 
936a7e14dcfSSatish Balay   ierr = MatZeroEntries(ipmP->K);CHKERRQ(ierr);
937a7e14dcfSSatish Balay   /* Copy H */
938a7e14dcfSSatish Balay   for (i=hstart;i<hend;i++) {
939a7e14dcfSSatish Balay     ierr = MatGetRow(tao->hessian,i,&ncols,&cols,&vals);CHKERRQ(ierr);
940a7e14dcfSSatish Balay     if (ncols > 0) {
941a7e14dcfSSatish Balay       ierr = MatSetValues(ipmP->K,1,&i,ncols,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
942a7e14dcfSSatish Balay     }
943a7e14dcfSSatish Balay     ierr = MatRestoreRow(tao->hessian,i,&ncols,&cols,&vals);CHKERRQ(ierr);
944a7e14dcfSSatish Balay   }
945a7e14dcfSSatish Balay 
946a7e14dcfSSatish Balay   /* Copy Ae and Ae' */
947a7e14dcfSSatish Balay   if (ipmP->me > 0) {
948a7e14dcfSSatish Balay     ierr = MatGetOwnershipRange(tao->jacobian_equality,&aestart,&aeend);CHKERRQ(ierr);
949a7e14dcfSSatish Balay     for (i=aestart;i<aeend;i++) {
950a7e14dcfSSatish Balay       ierr = MatGetRow(tao->jacobian_equality,i,&ncols,&cols,&vals);CHKERRQ(ierr);
951a7e14dcfSSatish Balay       if (ncols > 0) {
952a7e14dcfSSatish Balay         /*Ae*/
953a7e14dcfSSatish Balay         row = i+r1;
954a7e14dcfSSatish Balay         ierr = MatSetValues(ipmP->K,1,&row,ncols,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
955a7e14dcfSSatish Balay         /*Ae'*/
956a7e14dcfSSatish Balay         for (j=0;j<ncols;j++) {
957a7e14dcfSSatish Balay           newcol = i + c2;
958a7e14dcfSSatish Balay           newrow = cols[j];
959a7e14dcfSSatish Balay           newval = vals[j];
960a7e14dcfSSatish Balay           ierr = MatSetValues(ipmP->K,1,&newrow,1,&newcol,&newval,INSERT_VALUES);CHKERRQ(ierr);
961a7e14dcfSSatish Balay         }
962a7e14dcfSSatish Balay       }
963a7e14dcfSSatish Balay       ierr = MatRestoreRow(tao->jacobian_equality,i,&ncols,&cols,&vals);CHKERRQ(ierr);
964a7e14dcfSSatish Balay     }
965a7e14dcfSSatish Balay   }
966a7e14dcfSSatish Balay 
967a7e14dcfSSatish Balay   if (ipmP->nb > 0) {
968a7e14dcfSSatish Balay     ierr = MatGetOwnershipRange(ipmP->Ai,&aistart,&aiend);CHKERRQ(ierr);
969a7e14dcfSSatish Balay     /* Copy Ai,and Ai' */
970a7e14dcfSSatish Balay     for (i=aistart;i<aiend;i++) {
971a7e14dcfSSatish Balay       row = i+r2;
972a7e14dcfSSatish Balay       ierr = MatGetRow(ipmP->Ai,i,&ncols,&cols,&vals);CHKERRQ(ierr);
973a7e14dcfSSatish Balay       if (ncols > 0) {
974a7e14dcfSSatish Balay         /*Ai*/
975a7e14dcfSSatish Balay         ierr = MatSetValues(ipmP->K,1,&row,ncols,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
976a7e14dcfSSatish Balay         /*-Ai'*/
977a7e14dcfSSatish Balay         for (j=0;j<ncols;j++) {
978a7e14dcfSSatish Balay           newcol = i + c3;
979a7e14dcfSSatish Balay           newrow = cols[j];
980a7e14dcfSSatish Balay           newval = -vals[j];
981a7e14dcfSSatish Balay           ierr = MatSetValues(ipmP->K,1,&newrow,1,&newcol,&newval,INSERT_VALUES);CHKERRQ(ierr);
982a7e14dcfSSatish Balay         }
983a7e14dcfSSatish Balay       }
984a7e14dcfSSatish Balay       ierr = MatRestoreRow(ipmP->Ai,i,&ncols,&cols,&vals);CHKERRQ(ierr);
985a7e14dcfSSatish Balay     }
986a7e14dcfSSatish Balay 
987a7e14dcfSSatish Balay     /* -I */
988a7e14dcfSSatish Balay     for (i=kstart;i<kend;i++) {
989a7e14dcfSSatish Balay       if (i>=r2 && i<r3) {
990a7e14dcfSSatish Balay         newrow = i;
991a7e14dcfSSatish Balay         newcol = i-r2+c1;
992a7e14dcfSSatish Balay         newval = -1.0;
993a7e14dcfSSatish Balay         MatSetValues(ipmP->K,1,&newrow,1,&newcol,&newval,INSERT_VALUES);CHKERRQ(ierr);
994a7e14dcfSSatish Balay       }
995a7e14dcfSSatish Balay     }
996a7e14dcfSSatish Balay 
997a7e14dcfSSatish Balay     /* Copy L,Y */
998a7e14dcfSSatish Balay     ierr = VecGetOwnershipRange(ipmP->s,&sstart,&send);CHKERRQ(ierr);
9995e081366SBarry Smith     ierr = VecGetArrayRead(ipmP->lamdai,&l);CHKERRQ(ierr);
10005e081366SBarry Smith     ierr = VecGetArrayRead(ipmP->s,&y);CHKERRQ(ierr);
1001a7e14dcfSSatish Balay 
1002a7e14dcfSSatish Balay     for (i=sstart;i<send;i++) {
1003a7e14dcfSSatish Balay       newcols[0] = c1+i;
1004a7e14dcfSSatish Balay       newcols[1] = c3+i;
1005a7e14dcfSSatish Balay       newvals[0] = l[i-sstart];
1006a7e14dcfSSatish Balay       newvals[1] = y[i-sstart];
1007a7e14dcfSSatish Balay       newrow = r3+i;
1008a7e14dcfSSatish Balay       ierr = MatSetValues(ipmP->K,1,&newrow,2,newcols,newvals,INSERT_VALUES);CHKERRQ(ierr);
1009a7e14dcfSSatish Balay     }
1010a7e14dcfSSatish Balay 
10115e081366SBarry Smith     ierr = VecRestoreArrayRead(ipmP->lamdai,&l);CHKERRQ(ierr);
10125e081366SBarry Smith     ierr = VecRestoreArrayRead(ipmP->s,&y);CHKERRQ(ierr);
1013a7e14dcfSSatish Balay   }
1014a7e14dcfSSatish Balay 
1015a7e14dcfSSatish Balay   ierr = PetscFree(indices);CHKERRQ(ierr);
1016a7e14dcfSSatish Balay   ierr = PetscFree(newvals);CHKERRQ(ierr);
1017a7e14dcfSSatish Balay   ierr = MatAssemblyBegin(ipmP->K,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1018a7e14dcfSSatish Balay   ierr = MatAssemblyEnd(ipmP->K,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1019a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1020a7e14dcfSSatish Balay }
1021a7e14dcfSSatish Balay 
1022a7e14dcfSSatish Balay #undef __FUNCT__
1023a7e14dcfSSatish Balay #define __FUNCT__ "IPMGatherRHS"
1024441846f8SBarry Smith PetscErrorCode IPMGatherRHS(Tao tao,Vec RHS,Vec X1,Vec X2,Vec X3,Vec X4)
1025a7e14dcfSSatish Balay {
1026a7e14dcfSSatish Balay   TAO_IPM        *ipmP = (TAO_IPM *)tao->data;
1027a7e14dcfSSatish Balay   PetscErrorCode ierr;
1028a7e14dcfSSatish Balay 
102947a47007SBarry Smith   PetscFunctionBegin;
1030a7e14dcfSSatish Balay   /* rhs = [x1      (n)
1031a7e14dcfSSatish Balay             x2     (me)
1032a7e14dcfSSatish Balay             x3     (nb)
1033a7e14dcfSSatish Balay             x4     (nb)] */
1034a7e14dcfSSatish Balay   if (X1) {
1035a7e14dcfSSatish Balay     ierr = VecScatterBegin(ipmP->rhs1,X1,RHS,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1036a7e14dcfSSatish Balay     ierr = VecScatterEnd(ipmP->rhs1,X1,RHS,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1037a7e14dcfSSatish Balay   }
1038a7e14dcfSSatish Balay   if (ipmP->me > 0 && X2) {
1039a7e14dcfSSatish Balay     ierr = VecScatterBegin(ipmP->rhs2,X2,RHS,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1040a7e14dcfSSatish Balay     ierr = VecScatterEnd(ipmP->rhs2,X2,RHS,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1041a7e14dcfSSatish Balay   }
1042a7e14dcfSSatish Balay   if (ipmP->nb > 0) {
1043a7e14dcfSSatish Balay     if (X3) {
1044a7e14dcfSSatish Balay       ierr = VecScatterBegin(ipmP->rhs3,X3,RHS,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1045a7e14dcfSSatish Balay       ierr = VecScatterEnd(ipmP->rhs3,X3,RHS,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1046a7e14dcfSSatish Balay     }
1047a7e14dcfSSatish Balay     if (X4) {
1048a7e14dcfSSatish Balay       ierr = VecScatterBegin(ipmP->rhs4,X4,RHS,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1049a7e14dcfSSatish Balay       ierr = VecScatterEnd(ipmP->rhs4,X4,RHS,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1050a7e14dcfSSatish Balay     }
1051a7e14dcfSSatish Balay   }
1052a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1053a7e14dcfSSatish Balay }
1054a7e14dcfSSatish Balay 
1055a7e14dcfSSatish Balay #undef __FUNCT__
1056a7e14dcfSSatish Balay #define __FUNCT__ "IPMScatterStep"
1057441846f8SBarry Smith PetscErrorCode IPMScatterStep(Tao tao, Vec STEP, Vec X1, Vec X2, Vec X3, Vec X4)
1058a7e14dcfSSatish Balay {
1059a7e14dcfSSatish Balay   TAO_IPM        *ipmP = (TAO_IPM *)tao->data;
1060a7e14dcfSSatish Balay   PetscErrorCode ierr;
106147a47007SBarry Smith 
1062a7e14dcfSSatish Balay   PetscFunctionBegin;
1063a7e14dcfSSatish Balay   CHKMEMQ;
1064a7e14dcfSSatish Balay   /*        [x1    (n)
1065a7e14dcfSSatish Balay              x2    (nb) may be 0
1066a7e14dcfSSatish Balay              x3    (me) may be 0
1067a7e14dcfSSatish Balay              x4    (nb) may be 0 */
1068a7e14dcfSSatish Balay   if (X1) {
1069a7e14dcfSSatish Balay     ierr = VecScatterBegin(ipmP->step1,STEP,X1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1070a7e14dcfSSatish Balay     ierr = VecScatterEnd(ipmP->step1,STEP,X1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1071a7e14dcfSSatish Balay   }
1072a7e14dcfSSatish Balay   if (X2 && ipmP->nb > 0) {
1073a7e14dcfSSatish Balay     ierr = VecScatterBegin(ipmP->step2,STEP,X2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1074a7e14dcfSSatish Balay     ierr = VecScatterEnd(ipmP->step2,STEP,X2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1075a7e14dcfSSatish Balay   }
1076a7e14dcfSSatish Balay   if (X3 && ipmP->me > 0) {
1077a7e14dcfSSatish Balay     ierr = VecScatterBegin(ipmP->step3,STEP,X3,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1078a7e14dcfSSatish Balay     ierr = VecScatterEnd(ipmP->step3,STEP,X3,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1079a7e14dcfSSatish Balay   }
1080a7e14dcfSSatish Balay   if (X4 && ipmP->nb > 0) {
1081a7e14dcfSSatish Balay     ierr = VecScatterBegin(ipmP->step4,STEP,X4,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1082a7e14dcfSSatish Balay     ierr = VecScatterEnd(ipmP->step4,STEP,X4,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1083a7e14dcfSSatish Balay   }
1084a7e14dcfSSatish Balay   CHKMEMQ;
1085a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1086a7e14dcfSSatish Balay }
1087a7e14dcfSSatish Balay 
10881522df2eSJason Sarich /*MC
10891522df2eSJason Sarich   TAOIPM - Interior point algorithm for generally constrained optimization.
10901522df2eSJason Sarich 
10911522df2eSJason Sarich   Option Database Keys:
10921522df2eSJason Sarich +   -tao_ipm_pushnu - parameter to push initial dual variables away from bounds
10931522df2eSJason Sarich .   -tao_ipm_pushs - parameter to push initial slack variables away from bounds
10941522df2eSJason Sarich 
10951522df2eSJason Sarich   Notes: This algorithm is more of a place-holder for future constrained optimization algorithms and should not yet be used for large problems or production code.
10961eb8069cSJason Sarich   Level: beginner
10971eb8069cSJason Sarich 
10981522df2eSJason Sarich M*/
10991522df2eSJason Sarich 
1100a7e14dcfSSatish Balay #undef __FUNCT__
1101a7e14dcfSSatish Balay #define __FUNCT__ "TaoCreate_IPM"
1102728e0ed0SBarry Smith PETSC_EXTERN PetscErrorCode TaoCreate_IPM(Tao tao)
1103a7e14dcfSSatish Balay {
1104a7e14dcfSSatish Balay   TAO_IPM        *ipmP;
1105a7e14dcfSSatish Balay   PetscErrorCode ierr;
1106a7e14dcfSSatish Balay 
1107a7e14dcfSSatish Balay   PetscFunctionBegin;
1108a7e14dcfSSatish Balay   tao->ops->setup = TaoSetup_IPM;
1109a7e14dcfSSatish Balay   tao->ops->solve = TaoSolve_IPM;
1110a7e14dcfSSatish Balay   tao->ops->view = TaoView_IPM;
1111a7e14dcfSSatish Balay   tao->ops->setfromoptions = TaoSetFromOptions_IPM;
1112a7e14dcfSSatish Balay   tao->ops->destroy = TaoDestroy_IPM;
1113e9f9aeaeSSatish Balay   /* tao->ops->computedual = TaoComputeDual_IPM; */
1114a7e14dcfSSatish Balay 
11153c9e27cfSGeoffrey Irving   ierr = PetscNewLog(tao,&ipmP);CHKERRQ(ierr);
1116a7e14dcfSSatish Balay   tao->data = (void*)ipmP;
1117a7e14dcfSSatish Balay   tao->max_it = 200;
1118a7e14dcfSSatish Balay   tao->max_funcs = 500;
1119a7e14dcfSSatish Balay   tao->fatol = 1e-4;
1120a7e14dcfSSatish Balay   tao->frtol = 1e-4;
1121a7e14dcfSSatish Balay   ipmP->dec = 10000; /* line search critera */
1122a7e14dcfSSatish Balay   ipmP->taumin = 0.995;
1123a7e14dcfSSatish Balay   ipmP->monitorkkt = PETSC_FALSE;
1124a7e14dcfSSatish Balay   ipmP->pushs = 100;
1125a7e14dcfSSatish Balay   ipmP->pushnu = 100;
1126a7e14dcfSSatish Balay   ierr = KSPCreate(((PetscObject)tao)->comm, &tao->ksp);CHKERRQ(ierr);
1127a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1128a7e14dcfSSatish Balay }
1129728e0ed0SBarry Smith 
1130a7e14dcfSSatish Balay 
1131