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