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