xref: /petsc/src/tao/constrained/impls/ipm/ipm.c (revision aaa7dc30da3270cff6cb10b1db605b2ca746f216)
1*aaa7dc30SBarry Smith #include <taolinesearch.h>
2*aaa7dc30SBarry Smith #include <../src/tao/constrained/impls/ipm/ipm.h> /*I "ipm.h" I*/
3e9f9aeaeSSatish Balay /*
4e9f9aeaeSSatish Balay #define DEBUG_IPM
5e9f9aeaeSSatish Balay #define DEBUG_K
6e9f9aeaeSSatish Balay #define DEBUG_SCATTER
7e9f9aeaeSSatish Balay #define DEBUG_KKT
8e9f9aeaeSSatish Balay */
9a7e14dcfSSatish Balay /*
10a7e14dcfSSatish Balay    x,d in R^n
11a7e14dcfSSatish Balay    f in R
12a7e14dcfSSatish Balay    nb = mi + nlb+nub
13a7e14dcfSSatish Balay    s in R^nb is slack vector CI(x) / x-XL / -x+XU
14a7e14dcfSSatish Balay    bin in R^mi (tao->constraints_inequality)
15a7e14dcfSSatish Balay    beq in R^me (tao->constraints_equality)
16a7e14dcfSSatish Balay    lamdai in R^nb (ipmP->lamdai)
17a7e14dcfSSatish Balay    lamdae in R^me (ipmP->lamdae)
18a7e14dcfSSatish Balay    Jeq in R^(me x n) (tao->jacobian_equality)
19a7e14dcfSSatish Balay    Jin in R^(mi x n) (tao->jacobian_inequality)
20a7e14dcfSSatish Balay    Ai in  R^(nb x n) (ipmP->Ai)
21a7e14dcfSSatish Balay    H in R^(n x n) (tao->hessian)
22a7e14dcfSSatish Balay    min f=(1/2)*x'*H*x + d'*x
23a7e14dcfSSatish Balay    s.t.  CE(x) == 0
24a7e14dcfSSatish Balay          CI(x) >= 0
25a7e14dcfSSatish Balay          x >= tao->XL
26a7e14dcfSSatish Balay          -x >= -tao->XU
27a7e14dcfSSatish Balay */
28a7e14dcfSSatish Balay 
29a7e14dcfSSatish Balay static PetscErrorCode IPMComputeKKT(TaoSolver tao);
30a7e14dcfSSatish Balay static PetscErrorCode IPMPushInitialPoint(TaoSolver tao);
31a7e14dcfSSatish Balay static PetscErrorCode IPMEvaluate(TaoSolver tao);
32a7e14dcfSSatish Balay static PetscErrorCode IPMUpdateK(TaoSolver tao);
33a7e14dcfSSatish Balay static PetscErrorCode IPMUpdateAi(TaoSolver tao);
34a7e14dcfSSatish Balay static PetscErrorCode IPMGatherRHS(TaoSolver tao,Vec,Vec,Vec,Vec,Vec);
35a7e14dcfSSatish Balay static PetscErrorCode IPMScatterStep(TaoSolver tao,Vec,Vec,Vec,Vec,Vec);
36a7e14dcfSSatish Balay static PetscErrorCode IPMInitializeBounds(TaoSolver tao);
37a7e14dcfSSatish Balay 
38a7e14dcfSSatish Balay #undef __FUNCT__
39a7e14dcfSSatish Balay #define __FUNCT__ "TaoSolve_IPM"
40a7e14dcfSSatish Balay static PetscErrorCode TaoSolve_IPM(TaoSolver tao)
41a7e14dcfSSatish Balay {
42a7e14dcfSSatish Balay   PetscErrorCode             ierr;
43a7e14dcfSSatish Balay   TAO_IPM                    *ipmP = (TAO_IPM*)tao->data;
44a7e14dcfSSatish Balay   TaoSolverTerminationReason reason = TAO_CONTINUE_ITERATING;
45a7e14dcfSSatish Balay   PetscInt                   iter = 0,its,i;
46a7e14dcfSSatish Balay   PetscScalar                stepsize=1.0;
47a7e14dcfSSatish Balay   PetscScalar                step_s,step_l,alpha,tau,sigma,phi_target;
48a7e14dcfSSatish Balay 
4947a47007SBarry Smith   PetscFunctionBegin;
50a7e14dcfSSatish Balay   /* Push initial point away from bounds */
51a7e14dcfSSatish Balay #if defined(DEBUG_IPM)
52a7e14dcfSSatish Balay   ierr = VecNorm(tao->solution,NORM_2,&tau);CHKERRQ(ierr);
5347a47007SBarry Smith   ierr = VecView(tao->solution,0);CHKERRQ(ierr);
54a7e14dcfSSatish Balay #endif
55a7e14dcfSSatish Balay   ierr = IPMInitializeBounds(tao);CHKERRQ(ierr);
56a7e14dcfSSatish Balay   ierr = IPMPushInitialPoint(tao);CHKERRQ(ierr);
57a7e14dcfSSatish Balay #if defined(DEBUG_IPM)
58a7e14dcfSSatish Balay   ierr = VecNorm(tao->solution,NORM_2,&tau);CHKERRQ(ierr);
5947a47007SBarry Smith   ierr = VecView(tao->solution,0);CHKERRQ(ierr);
60e9f9aeaeSSatish Balay   /*  PetscPrintf(PETSC_COMM_WORLD,"||x0|| = %g\n",tau); */
61a7e14dcfSSatish Balay #endif
62a7e14dcfSSatish Balay   ierr = VecCopy(tao->solution,ipmP->rhs_x);CHKERRQ(ierr);
63a7e14dcfSSatish Balay   ierr = IPMEvaluate(tao);CHKERRQ(ierr);
64a7e14dcfSSatish Balay   ierr = IPMComputeKKT(tao);CHKERRQ(ierr);
65a7e14dcfSSatish Balay   ierr = TaoMonitor(tao,iter++,ipmP->kkt_f,ipmP->phi,0.0,1.0,&reason);
66a7e14dcfSSatish Balay 
67a7e14dcfSSatish Balay   while (reason == TAO_CONTINUE_ITERATING) {
68a7e14dcfSSatish Balay     ierr = IPMUpdateK(tao);CHKERRQ(ierr);
69a7e14dcfSSatish Balay     /*
70a7e14dcfSSatish Balay        rhs.x    = -rd
71a7e14dcfSSatish Balay        rhs.lame = -rpe
72a7e14dcfSSatish Balay        rhs.lami = -rpi
73a7e14dcfSSatish Balay        rhs.com  = -com
74a7e14dcfSSatish Balay     */
75a7e14dcfSSatish Balay 
76a7e14dcfSSatish Balay     ierr = VecCopy(ipmP->rd,ipmP->rhs_x);CHKERRQ(ierr);
77a7e14dcfSSatish Balay     if (ipmP->me > 0) {
78a7e14dcfSSatish Balay       ierr = VecCopy(ipmP->rpe,ipmP->rhs_lamdae);CHKERRQ(ierr);
79a7e14dcfSSatish Balay     }
80a7e14dcfSSatish Balay     if (ipmP->nb > 0) {
81a7e14dcfSSatish Balay         ierr = VecCopy(ipmP->rpi,ipmP->rhs_lamdai);CHKERRQ(ierr);
82a7e14dcfSSatish Balay         ierr = VecCopy(ipmP->complementarity,ipmP->rhs_s);CHKERRQ(ierr);
83a7e14dcfSSatish Balay     }
84a7e14dcfSSatish Balay     ierr = IPMGatherRHS(tao,ipmP->bigrhs,ipmP->rhs_x,ipmP->rhs_lamdae,
85a7e14dcfSSatish Balay                         ipmP->rhs_lamdai,ipmP->rhs_s);CHKERRQ(ierr);
86a7e14dcfSSatish Balay     ierr = VecScale(ipmP->bigrhs,-1.0);CHKERRQ(ierr);
87a7e14dcfSSatish Balay 
88a7e14dcfSSatish Balay     /* solve K * step = rhs */
89a7e14dcfSSatish Balay     ierr = KSPSetOperators(tao->ksp,ipmP->K,ipmP->K,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
90a7e14dcfSSatish Balay     ierr = KSPSolve(tao->ksp,ipmP->bigrhs,ipmP->bigstep);CHKERRQ(ierr);
91a7e14dcfSSatish Balay 
92a7e14dcfSSatish Balay     ierr = IPMScatterStep(tao,ipmP->bigstep,tao->stepdirection,ipmP->ds,
93a7e14dcfSSatish Balay                           ipmP->dlamdae,ipmP->dlamdai);CHKERRQ(ierr);
94a7e14dcfSSatish Balay     ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr);
95a7e14dcfSSatish Balay     tao->ksp_its += its;
96a7e14dcfSSatish Balay #if defined DEBUG_KKT
97a7e14dcfSSatish Balay     PetscPrintf(PETSC_COMM_WORLD,"first solve.\n");
98a7e14dcfSSatish Balay     PetscPrintf(PETSC_COMM_WORLD,"rhs_lamdai\n");
99e9f9aeaeSSatish Balay     /* VecView(ipmP->rhs_lamdai,0);
100e9f9aeaeSSatish Balay     ierr = VecView(ipmP->bigrhs,0);
101e9f9aeaeSSatish Balay     ierr = VecView(ipmP->bigstep,0); */
102a7e14dcfSSatish Balay     PetscScalar norm1,norm2;
103a7e14dcfSSatish Balay     ierr = VecNorm(ipmP->bigrhs,NORM_2,&norm1);
104a7e14dcfSSatish Balay     ierr = VecNorm(ipmP->bigstep,NORM_2,&norm2);
105a7e14dcfSSatish Balay     PetscPrintf(PETSC_COMM_WORLD,"||rhs|| = %g\t ||step|| = %g\n",norm1,norm2);
106a7e14dcfSSatish Balay #endif
107a7e14dcfSSatish Balay      /* Find distance along step direction to closest bound */
108a7e14dcfSSatish Balay     if (ipmP->nb > 0) {
1096c23d075SBarry Smith       ierr = VecStepBoundInfo(ipmP->s,ipmP->Zero_nb,ipmP->Inf_nb,ipmP->ds,&step_s,NULL,NULL);CHKERRQ(ierr);
1106c23d075SBarry Smith       ierr = VecStepBoundInfo(ipmP->lamdai,ipmP->Zero_nb,ipmP->Inf_nb,ipmP->dlamdai,&step_l,NULL,NULL);CHKERRQ(ierr);
111a7e14dcfSSatish Balay       alpha = PetscMin(step_s,step_l);
112a7e14dcfSSatish Balay       alpha = PetscMin(alpha,1.0);
113a7e14dcfSSatish Balay       ipmP->alpha1 = alpha;
114a7e14dcfSSatish Balay     } else {
115a7e14dcfSSatish Balay       ipmP->alpha1 = alpha = 1.0;
116a7e14dcfSSatish Balay     }
117a7e14dcfSSatish Balay 
118a7e14dcfSSatish Balay 
119a7e14dcfSSatish Balay     /* x_aff = x + alpha*d */
120a7e14dcfSSatish Balay     ierr = VecCopy(tao->solution,ipmP->save_x);CHKERRQ(ierr);
121a7e14dcfSSatish Balay     if (ipmP->me > 0) {
122a7e14dcfSSatish Balay       ierr = VecCopy(ipmP->lamdae,ipmP->save_lamdae);CHKERRQ(ierr);
123a7e14dcfSSatish Balay     }
124a7e14dcfSSatish Balay     if (ipmP->nb > 0) {
125a7e14dcfSSatish Balay       ierr = VecCopy(ipmP->lamdai,ipmP->save_lamdai);CHKERRQ(ierr);
126a7e14dcfSSatish Balay       ierr = VecCopy(ipmP->s,ipmP->save_s);CHKERRQ(ierr);
127a7e14dcfSSatish Balay     }
128a7e14dcfSSatish Balay 
129a7e14dcfSSatish Balay     ierr = VecAXPY(tao->solution,alpha,tao->stepdirection);CHKERRQ(ierr);
130a7e14dcfSSatish Balay     if (ipmP->me > 0) {
131a7e14dcfSSatish Balay       ierr = VecAXPY(ipmP->lamdae,alpha,ipmP->dlamdae);CHKERRQ(ierr);
132a7e14dcfSSatish Balay     }
133a7e14dcfSSatish Balay     if (ipmP->nb > 0) {
134a7e14dcfSSatish Balay       ierr = VecAXPY(ipmP->lamdai,alpha,ipmP->dlamdai);CHKERRQ(ierr);
135a7e14dcfSSatish Balay       ierr = VecAXPY(ipmP->s,alpha,ipmP->ds);CHKERRQ(ierr);
136a7e14dcfSSatish Balay     }
137a7e14dcfSSatish Balay 
138a7e14dcfSSatish Balay 
139a7e14dcfSSatish Balay     /* Recompute kkt to find centering parameter sigma = (new_mu/old_mu)^3 */
140a7e14dcfSSatish Balay     if (ipmP->mu == 0.0) {
141a7e14dcfSSatish Balay       sigma = 0.0;
142a7e14dcfSSatish Balay     } else {
143a7e14dcfSSatish Balay       sigma = 1.0/ipmP->mu;
144a7e14dcfSSatish Balay     }
145a7e14dcfSSatish Balay     ierr = IPMComputeKKT(tao);CHKERRQ(ierr);
146a7e14dcfSSatish Balay     sigma *= ipmP->mu;
147a7e14dcfSSatish Balay     sigma*=sigma*sigma;
148a7e14dcfSSatish Balay 
149a7e14dcfSSatish Balay     /* revert kkt info */
150a7e14dcfSSatish Balay     ierr = VecCopy(ipmP->save_x,tao->solution);CHKERRQ(ierr);
151a7e14dcfSSatish Balay     if (ipmP->me > 0) {
152a7e14dcfSSatish Balay       ierr = VecCopy(ipmP->save_lamdae,ipmP->lamdae);CHKERRQ(ierr);
153a7e14dcfSSatish Balay     }
154a7e14dcfSSatish Balay     if (ipmP->nb > 0) {
155a7e14dcfSSatish Balay       ierr = VecCopy(ipmP->save_lamdai,ipmP->lamdai);CHKERRQ(ierr);
156a7e14dcfSSatish Balay       ierr = VecCopy(ipmP->save_s,ipmP->s);CHKERRQ(ierr);
157a7e14dcfSSatish Balay     }
158a7e14dcfSSatish Balay     ierr = IPMComputeKKT(tao);CHKERRQ(ierr);
159a7e14dcfSSatish Balay 
160a7e14dcfSSatish Balay     /* update rhs with new complementarity vector */
161a7e14dcfSSatish Balay     if (ipmP->nb > 0) {
162a7e14dcfSSatish Balay       ierr = VecCopy(ipmP->complementarity,ipmP->rhs_s);CHKERRQ(ierr);
163a7e14dcfSSatish Balay       ierr = VecScale(ipmP->rhs_s,-1.0);CHKERRQ(ierr);
164a7e14dcfSSatish Balay       ierr = VecShift(ipmP->rhs_s,sigma*ipmP->mu);CHKERRQ(ierr);
165a7e14dcfSSatish Balay     }
1666c23d075SBarry Smith     ierr = IPMGatherRHS(tao,ipmP->bigrhs,NULL,NULL,NULL,ipmP->rhs_s);CHKERRQ(ierr);
167a7e14dcfSSatish Balay 
168a7e14dcfSSatish Balay     /* solve K * step = rhs */
169a7e14dcfSSatish Balay     ierr = KSPSetOperators(tao->ksp,ipmP->K,ipmP->K,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
170a7e14dcfSSatish Balay     ierr = KSPSolve(tao->ksp,ipmP->bigrhs,ipmP->bigstep);CHKERRQ(ierr);
171a7e14dcfSSatish Balay 
17247a47007SBarry Smith     ierr = IPMScatterStep(tao,ipmP->bigstep,tao->stepdirection,ipmP->ds,ipmP->dlamdae,ipmP->dlamdai);CHKERRQ(ierr);
173a7e14dcfSSatish Balay     ierr = KSPGetIterationNumber(tao->ksp,&its);CHKERRQ(ierr);
174a7e14dcfSSatish Balay #if defined DEBUG_KKT2
175a7e14dcfSSatish Balay     PetscPrintf(PETSC_COMM_WORLD,"rhs_lamdai\n");
176a7e14dcfSSatish Balay     VecView(ipmP->rhs_lamdai,0);
177a7e14dcfSSatish Balay     ierr = VecView(ipmP->bigrhs,0);
178a7e14dcfSSatish Balay     ierr = VecView(ipmP->bigstep,0);
179a7e14dcfSSatish Balay #endif
180a7e14dcfSSatish Balay     tao->ksp_its += its;
181a7e14dcfSSatish Balay 
182a7e14dcfSSatish Balay     if (ipmP->nb > 0) {
183a7e14dcfSSatish Balay       /* Get max step size and apply frac-to-boundary */
184a7e14dcfSSatish Balay       tau = PetscMax(ipmP->taumin,1.0-ipmP->mu);
185a7e14dcfSSatish Balay       tau = PetscMin(tau,1.0);
186a7e14dcfSSatish Balay       if (tau != 1.0) {
187a7e14dcfSSatish Balay         ierr = VecScale(ipmP->s,tau);CHKERRQ(ierr);
188a7e14dcfSSatish Balay         ierr = VecScale(ipmP->lamdai,tau);CHKERRQ(ierr);
189a7e14dcfSSatish Balay       }
1906c23d075SBarry Smith       ierr = VecStepBoundInfo(ipmP->s,ipmP->Zero_nb,ipmP->Inf_nb,ipmP->ds,&step_s,NULL,NULL);CHKERRQ(ierr);
1916c23d075SBarry Smith       ierr = VecStepBoundInfo(ipmP->lamdai,ipmP->Zero_nb,ipmP->Inf_nb,ipmP->dlamdai,&step_l,NULL,NULL);CHKERRQ(ierr);
192a7e14dcfSSatish Balay       if (tau != 1.0) {
193a7e14dcfSSatish Balay         ierr = VecCopy(ipmP->save_s,ipmP->s);CHKERRQ(ierr);
194a7e14dcfSSatish Balay         ierr = VecCopy(ipmP->save_lamdai,ipmP->lamdai);CHKERRQ(ierr);
195a7e14dcfSSatish Balay       }
196a7e14dcfSSatish Balay       alpha = PetscMin(step_s,step_l);
197a7e14dcfSSatish Balay       alpha = PetscMin(alpha,1.0);
198a7e14dcfSSatish Balay     } else {
199a7e14dcfSSatish Balay       alpha = 1.0;
200a7e14dcfSSatish Balay     }
201a7e14dcfSSatish Balay     ipmP->alpha2 = alpha;
202a7e14dcfSSatish Balay     /* TODO make phi_target meaningful */
203a7e14dcfSSatish Balay     phi_target = ipmP->dec * ipmP->phi;
204a7e14dcfSSatish Balay     for (i=0; i<11;i++) {
205a7e14dcfSSatish Balay #if defined DEBUG_KKT2
206a7e14dcfSSatish Balay     PetscPrintf(PETSC_COMM_WORLD,"alpha2=%g\n",alpha);
207a7e14dcfSSatish Balay       PetscPrintf(PETSC_COMM_WORLD,"old point:\n");
208a7e14dcfSSatish Balay       VecView(tao->solution,0);
209a7e14dcfSSatish Balay       VecView(ipmP->lamdae,0);
210a7e14dcfSSatish Balay       VecView(ipmP->s,0);
211a7e14dcfSSatish Balay       VecView(ipmP->lamdai,0);
212a7e14dcfSSatish Balay #endif
213a7e14dcfSSatish Balay       ierr = VecAXPY(tao->solution,alpha,tao->stepdirection);CHKERRQ(ierr);
214a7e14dcfSSatish Balay       if (ipmP->nb > 0) {
215a7e14dcfSSatish Balay         ierr = VecAXPY(ipmP->s,alpha,ipmP->ds);CHKERRQ(ierr);
216a7e14dcfSSatish Balay         ierr = VecAXPY(ipmP->lamdai,alpha,ipmP->dlamdai);CHKERRQ(ierr);
217a7e14dcfSSatish Balay       }
218a7e14dcfSSatish Balay       if (ipmP->me > 0) {
219a7e14dcfSSatish Balay         ierr = VecAXPY(ipmP->lamdae,alpha,ipmP->dlamdae);CHKERRQ(ierr);
220a7e14dcfSSatish Balay       }
221a7e14dcfSSatish Balay #if defined DEBUG_KKT
222a7e14dcfSSatish Balay       PetscPrintf(PETSC_COMM_WORLD,"step direction:\n");
223a7e14dcfSSatish Balay       VecView(tao->stepdirection,0);
224e9f9aeaeSSatish Balay       /* VecView(ipmP->dlamdae,0);
225e9f9aeaeSSatish Balay       VecView(ipmP->ds,0);
226e9f9aeaeSSatish Balay       VecView(ipmP->dlamdai,0);
227a7e14dcfSSatish Balay 
228e9f9aeaeSSatish Balay       PetscPrintf(PETSC_COMM_WORLD,"New iterate:\n");
229e9f9aeaeSSatish Balay       VecView(tao->solution,0);
230e9f9aeaeSSatish Balay       VecView(ipmP->lamdae,0);
231e9f9aeaeSSatish Balay       VecView(ipmP->s,0);
232e9f9aeaeSSatish Balay       VecView(ipmP->lamdai,0); */
233a7e14dcfSSatish Balay #endif
234a7e14dcfSSatish Balay       /* update dual variables */
235a7e14dcfSSatish Balay       if (ipmP->me > 0) {
236a7e14dcfSSatish Balay         ierr = VecCopy(ipmP->lamdae,tao->DE);CHKERRQ(ierr);
237a7e14dcfSSatish Balay       }
238a7e14dcfSSatish Balay       /* TODO: fix
239a7e14dcfSSatish Balay       if (ipmP->nb > 0) {
240a7e14dcfSSatish Balay         ierr = VecScatterBegin
241a7e14dcfSSatish Balay         PetscInt lstart,lend;
242a7e14dcfSSatish Balay 
243a7e14dcfSSatish Balay         ierr = VecGetOwnershipRange(ipmP->lamdai,&lstart,&lend);
244a7e14dcfSSatish Balay         ierr = VecGetArray(ipmP->lamdai,&li);CHKERRQ(ierr);
245a7e14dcfSSatish Balay         ierr = VecGetArray(tao->DI,&di);CHKERRQ(ierr);
246a7e14dcfSSatish Balay         for (j=lstart;j<lend;j++) {
247a7e14dcfSSatish Balay           if (j < ipmP->nilb) {
248a7e14dcfSSatish Balay             di[j] = li[j];
249a7e14dcfSSatish Balay           }
250a7e14dcfSSatish Balay         }
251a7e14dcfSSatish Balay         ierr = VecRestoreArray(ipmP->lamdai,&li);CHKERRQ(ierr);
252a7e14dcfSSatish Balay         ierr = VecRestoreArray(tao->DI,&di);CHKERRQ(ierr);
253a7e14dcfSSatish Balay       }
254a7e14dcfSSatish Balay       */
255a7e14dcfSSatish Balay 
256a7e14dcfSSatish Balay 
257a7e14dcfSSatish Balay       ierr = IPMEvaluate(tao);CHKERRQ(ierr);
258a7e14dcfSSatish Balay       ierr = IPMComputeKKT(tao);CHKERRQ(ierr);
259a7e14dcfSSatish Balay       if (ipmP->phi <= phi_target) break;
260a7e14dcfSSatish Balay       alpha /= 2.0;
261a7e14dcfSSatish Balay     }
262a7e14dcfSSatish Balay 
26347a47007SBarry Smith     ierr = TaoMonitor(tao,iter,ipmP->kkt_f,ipmP->phi,0.0,stepsize,&reason);CHKERRQ(ierr);
264a7e14dcfSSatish Balay     iter++;
265a7e14dcfSSatish Balay   }
266a7e14dcfSSatish Balay   PetscFunctionReturn(0);
267a7e14dcfSSatish Balay }
268a7e14dcfSSatish Balay 
269a7e14dcfSSatish Balay #undef __FUNCT__
270a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetup_IPM"
271a7e14dcfSSatish Balay static PetscErrorCode TaoSetup_IPM(TaoSolver tao)
272a7e14dcfSSatish Balay {
273a7e14dcfSSatish Balay   TAO_IPM        *ipmP = (TAO_IPM*)tao->data;
274a7e14dcfSSatish Balay   PetscErrorCode ierr;
275a7e14dcfSSatish Balay 
276a7e14dcfSSatish Balay   PetscFunctionBegin;
277a7e14dcfSSatish Balay   ipmP->nb = ipmP->mi = ipmP->me = 0;
278a7e14dcfSSatish Balay   ipmP->K=0;
279a7e14dcfSSatish Balay   ierr = VecGetSize(tao->solution,&ipmP->n);CHKERRQ(ierr);
280a7e14dcfSSatish Balay   if (!tao->gradient) {
281a7e14dcfSSatish Balay     ierr = VecDuplicate(tao->solution, &tao->gradient);CHKERRQ(ierr);
282a7e14dcfSSatish Balay     ierr = VecDuplicate(tao->solution, &tao->stepdirection);CHKERRQ(ierr);
283a7e14dcfSSatish Balay     ierr = VecDuplicate(tao->solution, &ipmP->rd);CHKERRQ(ierr);
284a7e14dcfSSatish Balay     ierr = VecDuplicate(tao->solution, &ipmP->rhs_x);CHKERRQ(ierr);
285a7e14dcfSSatish Balay     ierr = VecDuplicate(tao->solution, &ipmP->work);CHKERRQ(ierr);
286a7e14dcfSSatish Balay     ierr = VecDuplicate(tao->solution, &ipmP->save_x);CHKERRQ(ierr);
287a7e14dcfSSatish Balay   }
288a7e14dcfSSatish Balay   if (tao->constraints_equality) {
289a7e14dcfSSatish Balay     ierr = VecGetSize(tao->constraints_equality,&ipmP->me);CHKERRQ(ierr);
290a7e14dcfSSatish Balay     ierr = VecDuplicate(tao->constraints_equality,&ipmP->lamdae);CHKERRQ(ierr);
291a7e14dcfSSatish Balay     ierr = VecDuplicate(tao->constraints_equality,&ipmP->dlamdae);CHKERRQ(ierr);
292a7e14dcfSSatish Balay     ierr = VecDuplicate(tao->constraints_equality,&ipmP->rhs_lamdae);CHKERRQ(ierr);
293a7e14dcfSSatish Balay     ierr = VecDuplicate(tao->constraints_equality,&ipmP->save_lamdae);CHKERRQ(ierr);
294a7e14dcfSSatish Balay     ierr = VecDuplicate(tao->constraints_equality,&ipmP->rpe);CHKERRQ(ierr);
295a7e14dcfSSatish Balay     ierr = VecDuplicate(tao->constraints_equality,&tao->DE);CHKERRQ(ierr);
296a7e14dcfSSatish Balay   }
297a7e14dcfSSatish Balay   if (tao->constraints_inequality) {
298a7e14dcfSSatish Balay     ierr = VecDuplicate(tao->constraints_inequality,&tao->DI);CHKERRQ(ierr);
299a7e14dcfSSatish Balay   }
300a7e14dcfSSatish Balay   PetscFunctionReturn(0);
301a7e14dcfSSatish Balay }
302a7e14dcfSSatish Balay 
303a7e14dcfSSatish Balay #undef __FUNCT__
304a7e14dcfSSatish Balay #define __FUNCT__ "IPMInitializeBounds"
305a7e14dcfSSatish Balay static PetscErrorCode IPMInitializeBounds(TaoSolver tao)
306a7e14dcfSSatish Balay {
307a7e14dcfSSatish Balay   TAO_IPM        *ipmP = (TAO_IPM*)tao->data;
308a7e14dcfSSatish Balay   Vec            xtmp;
309a7e14dcfSSatish Balay   PetscInt       xstart,xend;
310a7e14dcfSSatish Balay   PetscInt       ucstart,ucend; /* user ci */
311a7e14dcfSSatish Balay   PetscInt       ucestart,uceend; /* user ce */
312a7e14dcfSSatish Balay   PetscInt       sstart,send;
313a7e14dcfSSatish Balay   PetscInt       bigsize;
314a7e14dcfSSatish Balay   PetscInt       i,counter,nloc;
315a7e14dcfSSatish Balay   PetscInt       *cind,*xind,*ucind,*uceind,*stepind;
316a7e14dcfSSatish Balay   VecType        vtype;
317a7e14dcfSSatish Balay   const PetscInt *xli,*xui;
318a7e14dcfSSatish Balay   PetscInt       xl_offset,xu_offset;
319a7e14dcfSSatish Balay   IS             bigxl,bigxu,isuc,isc,isx,sis,is1;
320a7e14dcfSSatish Balay   PetscErrorCode ierr;
321a7e14dcfSSatish Balay   MPI_Comm       comm;
32247a47007SBarry Smith 
323a7e14dcfSSatish Balay   PetscFunctionBegin;
324a7e14dcfSSatish Balay   cind=xind=ucind=uceind=stepind=0;
325a7e14dcfSSatish Balay   ipmP->mi=0;
326a7e14dcfSSatish Balay   ipmP->nxlb=0;
327a7e14dcfSSatish Balay   ipmP->nxub=0;
328a7e14dcfSSatish Balay   ipmP->nb=0;
329a7e14dcfSSatish Balay   ipmP->nslack=0;
330a7e14dcfSSatish Balay 
331a7e14dcfSSatish Balay   ierr = VecDuplicate(tao->solution,&xtmp);CHKERRQ(ierr);
332a7e14dcfSSatish Balay   if (!tao->XL && !tao->XU && tao->ops->computebounds) {
333a7e14dcfSSatish Balay     ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr);
334a7e14dcfSSatish Balay   }
335a7e14dcfSSatish Balay   if (tao->XL) {
336a7e14dcfSSatish Balay     ierr = VecSet(xtmp,TAO_NINFINITY);CHKERRQ(ierr);
337a7e14dcfSSatish Balay     ierr = VecWhichGreaterThan(tao->XL,xtmp,&ipmP->isxl);CHKERRQ(ierr);
338a7e14dcfSSatish Balay     ierr = ISGetSize(ipmP->isxl,&ipmP->nxlb);CHKERRQ(ierr);
339a7e14dcfSSatish Balay   } else {
340a7e14dcfSSatish Balay     ipmP->nxlb=0;
341a7e14dcfSSatish Balay   }
342a7e14dcfSSatish Balay   if (tao->XU) {
343a7e14dcfSSatish Balay     ierr = VecSet(xtmp,TAO_INFINITY);CHKERRQ(ierr);
344a7e14dcfSSatish Balay     ierr = VecWhichLessThan(tao->XU,xtmp,&ipmP->isxu);CHKERRQ(ierr);
345a7e14dcfSSatish Balay     ierr = ISGetSize(ipmP->isxu,&ipmP->nxub);CHKERRQ(ierr);
346a7e14dcfSSatish Balay   } else {
347a7e14dcfSSatish Balay     ipmP->nxub=0;
348a7e14dcfSSatish Balay   }
349a7e14dcfSSatish Balay   ierr = VecDestroy(&xtmp);CHKERRQ(ierr);
350a7e14dcfSSatish Balay   if (tao->constraints_inequality) {
351a7e14dcfSSatish Balay     ierr = VecGetSize(tao->constraints_inequality,&ipmP->mi);CHKERRQ(ierr);
352a7e14dcfSSatish Balay   } else {
353a7e14dcfSSatish Balay     ipmP->mi = 0;
354a7e14dcfSSatish Balay   }
355a7e14dcfSSatish Balay #if defined DEBUG_K
356a7e14dcfSSatish Balay   PetscPrintf(PETSC_COMM_WORLD,"isxl:\n");
357a7e14dcfSSatish Balay   if (ipmP->nxlb) {
358a7e14dcfSSatish Balay     ISView(ipmP->isxl,0);
359a7e14dcfSSatish Balay   }
360a7e14dcfSSatish Balay   PetscPrintf(PETSC_COMM_WORLD,"isxu:\n");
361a7e14dcfSSatish Balay   if (ipmP->nxub) {
362a7e14dcfSSatish Balay     ISView(ipmP->isxu,0);
363a7e14dcfSSatish Balay   }
364a7e14dcfSSatish Balay 
365a7e14dcfSSatish Balay #endif
366a7e14dcfSSatish Balay   ipmP->nb = ipmP->nxlb + ipmP->nxub + ipmP->mi;
367a7e14dcfSSatish Balay 
368a7e14dcfSSatish Balay   comm = ((PetscObject)(tao->solution))->comm;
369a7e14dcfSSatish Balay 
370a7e14dcfSSatish Balay   bigsize = ipmP->n+2*ipmP->nb+ipmP->me;
3710e660641SBarry Smith   ierr = PetscMalloc1(bigsize,&stepind);CHKERRQ(ierr);
3720e660641SBarry Smith   ierr = PetscMalloc1(ipmP->n,&xind);CHKERRQ(ierr);
3730e660641SBarry Smith   ierr = PetscMalloc1(ipmP->me,&uceind);CHKERRQ(ierr);
374a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(tao->solution,&xstart,&xend);CHKERRQ(ierr);
375a7e14dcfSSatish Balay 
376a7e14dcfSSatish Balay   if (ipmP->nb > 0) {
377a7e14dcfSSatish Balay     ierr = VecCreate(comm,&ipmP->s);CHKERRQ(ierr);
378a7e14dcfSSatish Balay     ierr = VecSetSizes(ipmP->s,PETSC_DECIDE,ipmP->nb);CHKERRQ(ierr);
379a7e14dcfSSatish Balay     ierr = VecSetFromOptions(ipmP->s);CHKERRQ(ierr);
380a7e14dcfSSatish Balay     ierr = VecDuplicate(ipmP->s,&ipmP->ds);CHKERRQ(ierr);
381a7e14dcfSSatish Balay     ierr = VecDuplicate(ipmP->s,&ipmP->rhs_s);CHKERRQ(ierr);
382a7e14dcfSSatish Balay     ierr = VecDuplicate(ipmP->s,&ipmP->complementarity);CHKERRQ(ierr);
383a7e14dcfSSatish Balay     ierr = VecDuplicate(ipmP->s,&ipmP->ci);CHKERRQ(ierr);
384a7e14dcfSSatish Balay 
385a7e14dcfSSatish Balay     ierr = VecDuplicate(ipmP->s,&ipmP->lamdai);CHKERRQ(ierr);
386a7e14dcfSSatish Balay     ierr = VecDuplicate(ipmP->s,&ipmP->dlamdai);CHKERRQ(ierr);
387a7e14dcfSSatish Balay     ierr = VecDuplicate(ipmP->s,&ipmP->rhs_lamdai);CHKERRQ(ierr);
388a7e14dcfSSatish Balay     ierr = VecDuplicate(ipmP->s,&ipmP->save_lamdai);CHKERRQ(ierr);
389a7e14dcfSSatish Balay 
390a7e14dcfSSatish Balay     ierr = VecDuplicate(ipmP->s,&ipmP->save_s);CHKERRQ(ierr);
391a7e14dcfSSatish Balay     ierr = VecDuplicate(ipmP->s,&ipmP->rpi);CHKERRQ(ierr);
392a7e14dcfSSatish Balay     ierr = VecDuplicate(ipmP->s,&ipmP->Zero_nb);CHKERRQ(ierr);
393a7e14dcfSSatish Balay     ierr = VecSet(ipmP->Zero_nb,0.0);CHKERRQ(ierr);
394a7e14dcfSSatish Balay     ierr = VecDuplicate(ipmP->s,&ipmP->One_nb);CHKERRQ(ierr);
395a7e14dcfSSatish Balay     ierr = VecSet(ipmP->One_nb,1.0);CHKERRQ(ierr);
396a7e14dcfSSatish Balay     ierr = VecDuplicate(ipmP->s,&ipmP->Inf_nb);CHKERRQ(ierr);
397a7e14dcfSSatish Balay     ierr = VecSet(ipmP->Inf_nb,TAO_INFINITY);CHKERRQ(ierr);
398a7e14dcfSSatish Balay 
3990e660641SBarry Smith     ierr = PetscMalloc1(ipmP->nb,&cind);CHKERRQ(ierr);
4000e660641SBarry Smith     ierr = PetscMalloc1(ipmP->mi,&ucind);CHKERRQ(ierr);
401a7e14dcfSSatish Balay     ierr = VecGetOwnershipRange(ipmP->s,&sstart,&send);CHKERRQ(ierr);
402a7e14dcfSSatish Balay 
403a7e14dcfSSatish Balay     if (ipmP->mi > 0) {
404a7e14dcfSSatish Balay       ierr = VecGetOwnershipRange(tao->constraints_inequality,&ucstart,&ucend);CHKERRQ(ierr);
405a7e14dcfSSatish Balay       counter=0;
406a7e14dcfSSatish Balay       for (i=ucstart;i<ucend;i++) {
407a7e14dcfSSatish Balay         cind[counter++] = i;
408a7e14dcfSSatish Balay       }
409a7e14dcfSSatish Balay       ierr = ISCreateGeneral(comm,counter,cind,PETSC_COPY_VALUES,&isuc);CHKERRQ(ierr);
410a7e14dcfSSatish Balay       ierr = ISCreateGeneral(comm,counter,cind,PETSC_COPY_VALUES,&isc);CHKERRQ(ierr);
411a7e14dcfSSatish Balay       ierr = VecScatterCreate(tao->constraints_inequality,isuc,ipmP->ci,isc,&ipmP->ci_scat);CHKERRQ(ierr);
412a7e14dcfSSatish Balay 
413a7e14dcfSSatish Balay       ierr = ISDestroy(&isuc);
414a7e14dcfSSatish Balay       ierr = ISDestroy(&isc);
415a7e14dcfSSatish Balay     }
416a7e14dcfSSatish Balay     /* need to know how may xbound indices are on each process */
417a7e14dcfSSatish Balay     /* TODO better way */
418a7e14dcfSSatish Balay     if (ipmP->nxlb) {
419a7e14dcfSSatish Balay       ierr = ISAllGather(ipmP->isxl,&bigxl);CHKERRQ(ierr);
420a7e14dcfSSatish Balay       ierr = ISGetIndices(bigxl,&xli);CHKERRQ(ierr);
421a7e14dcfSSatish Balay       /* find offsets for this processor */
422a7e14dcfSSatish Balay       xl_offset = ipmP->mi;
423a7e14dcfSSatish Balay       for (i=0;i<ipmP->nxlb;i++) {
424a7e14dcfSSatish Balay         if (xli[i] < xstart) {
425a7e14dcfSSatish Balay           xl_offset++;
426a7e14dcfSSatish Balay         } else break;
427a7e14dcfSSatish Balay       }
428a7e14dcfSSatish Balay       ierr = ISRestoreIndices(bigxl,&xli);CHKERRQ(ierr);
429a7e14dcfSSatish Balay 
430a7e14dcfSSatish Balay       ierr = ISGetIndices(ipmP->isxl,&xli);CHKERRQ(ierr);
431a7e14dcfSSatish Balay       ierr = ISGetLocalSize(ipmP->isxl,&nloc);CHKERRQ(ierr);
432a7e14dcfSSatish Balay       for (i=0;i<nloc;i++) {
433a7e14dcfSSatish Balay         xind[i] = xli[i];
434a7e14dcfSSatish Balay         cind[i] = xl_offset+i;
435a7e14dcfSSatish Balay       }
436a7e14dcfSSatish Balay 
437a7e14dcfSSatish Balay       ierr = ISCreateGeneral(comm,nloc,xind,PETSC_COPY_VALUES,&isx);CHKERRQ(ierr);
438a7e14dcfSSatish Balay       ierr = ISCreateGeneral(comm,nloc,cind,PETSC_COPY_VALUES,&isc);CHKERRQ(ierr);
439a7e14dcfSSatish Balay       ierr = VecScatterCreate(tao->XL,isx,ipmP->ci,isc,&ipmP->xl_scat);CHKERRQ(ierr);
440a7e14dcfSSatish Balay       ierr = ISDestroy(&isx);CHKERRQ(ierr);
441a7e14dcfSSatish Balay       ierr = ISDestroy(&isc);CHKERRQ(ierr);
442a7e14dcfSSatish Balay       ierr = ISDestroy(&bigxl);CHKERRQ(ierr);
443a7e14dcfSSatish Balay     }
444a7e14dcfSSatish Balay 
445a7e14dcfSSatish Balay     if (ipmP->nxub) {
446a7e14dcfSSatish Balay       ierr = ISAllGather(ipmP->isxu,&bigxu);CHKERRQ(ierr);
447a7e14dcfSSatish Balay       ierr = ISGetIndices(bigxu,&xui);CHKERRQ(ierr);
448a7e14dcfSSatish Balay       /* find offsets for this processor */
449a7e14dcfSSatish Balay       xu_offset = ipmP->mi + ipmP->nxlb;
450a7e14dcfSSatish Balay       for (i=0;i<ipmP->nxub;i++) {
451a7e14dcfSSatish Balay         if (xui[i] < xstart) {
452a7e14dcfSSatish Balay           xu_offset++;
453a7e14dcfSSatish Balay         } else break;
454a7e14dcfSSatish Balay       }
455a7e14dcfSSatish Balay       ierr = ISRestoreIndices(bigxu,&xui);CHKERRQ(ierr);
456a7e14dcfSSatish Balay 
457a7e14dcfSSatish Balay       ierr = ISGetIndices(ipmP->isxu,&xui);CHKERRQ(ierr);
458a7e14dcfSSatish Balay       ierr = ISGetLocalSize(ipmP->isxu,&nloc);CHKERRQ(ierr);
459a7e14dcfSSatish Balay       for (i=0;i<nloc;i++) {
460a7e14dcfSSatish Balay         xind[i] = xui[i];
461a7e14dcfSSatish Balay         cind[i] = xu_offset+i;
462a7e14dcfSSatish Balay       }
463a7e14dcfSSatish Balay 
464a7e14dcfSSatish Balay       ierr = ISCreateGeneral(comm,nloc,xind,PETSC_COPY_VALUES,&isx);CHKERRQ(ierr);
465a7e14dcfSSatish Balay       ierr = ISCreateGeneral(comm,nloc,cind,PETSC_COPY_VALUES,&isc);CHKERRQ(ierr);
466a7e14dcfSSatish Balay       ierr = VecScatterCreate(tao->XU,isx,ipmP->ci,isc,&ipmP->xu_scat);CHKERRQ(ierr);
467a7e14dcfSSatish Balay       ierr = ISDestroy(&isx);CHKERRQ(ierr);
468a7e14dcfSSatish Balay       ierr = ISDestroy(&isc);CHKERRQ(ierr);
469a7e14dcfSSatish Balay       ierr = ISDestroy(&bigxu);CHKERRQ(ierr);
470a7e14dcfSSatish Balay     }
471a7e14dcfSSatish Balay   }
472a7e14dcfSSatish Balay   ierr = VecCreate(comm,&ipmP->bigrhs);CHKERRQ(ierr);
473a7e14dcfSSatish Balay   ierr = VecGetType(tao->solution,&vtype);CHKERRQ(ierr);
474a7e14dcfSSatish Balay   ierr = VecSetType(ipmP->bigrhs,vtype);CHKERRQ(ierr);
475a7e14dcfSSatish Balay   ierr = VecSetSizes(ipmP->bigrhs,PETSC_DECIDE,bigsize);CHKERRQ(ierr);
476a7e14dcfSSatish Balay   ierr = VecSetFromOptions(ipmP->bigrhs);CHKERRQ(ierr);
477a7e14dcfSSatish Balay   ierr = VecDuplicate(ipmP->bigrhs,&ipmP->bigstep);CHKERRQ(ierr);
478a7e14dcfSSatish Balay 
479a7e14dcfSSatish Balay   /* create scatters for step->x and x->rhs */
480a7e14dcfSSatish Balay   for (i=xstart;i<xend;i++) {
481a7e14dcfSSatish Balay     stepind[i-xstart] = i;
482a7e14dcfSSatish Balay     xind[i-xstart] = i;
483a7e14dcfSSatish Balay   }
484a7e14dcfSSatish Balay   ierr = ISCreateGeneral(comm,xend-xstart,stepind,PETSC_COPY_VALUES,&sis);CHKERRQ(ierr);
485a7e14dcfSSatish Balay   ierr = ISCreateGeneral(comm,xend-xstart,xind,PETSC_COPY_VALUES,&is1);CHKERRQ(ierr);
486a7e14dcfSSatish Balay   ierr = VecScatterCreate(ipmP->bigstep,sis,tao->solution,is1,&ipmP->step1);CHKERRQ(ierr);
487a7e14dcfSSatish Balay   ierr = VecScatterCreate(tao->solution,is1,ipmP->bigrhs,sis,&ipmP->rhs1);CHKERRQ(ierr);
488a7e14dcfSSatish Balay   ierr = ISDestroy(&sis);CHKERRQ(ierr);
489a7e14dcfSSatish Balay   ierr = ISDestroy(&is1);CHKERRQ(ierr);
490a7e14dcfSSatish Balay 
491a7e14dcfSSatish Balay   if (ipmP->nb > 0) {
492a7e14dcfSSatish Balay     for (i=sstart;i<send;i++) {
493a7e14dcfSSatish Balay       stepind[i-sstart] = i+ipmP->n;
494a7e14dcfSSatish Balay       cind[i-sstart] = i;
495a7e14dcfSSatish Balay     }
496a7e14dcfSSatish Balay     ierr = ISCreateGeneral(comm,send-sstart,stepind,PETSC_COPY_VALUES,&sis);CHKERRQ(ierr);
497a7e14dcfSSatish Balay     ierr = ISCreateGeneral(comm,send-sstart,cind,PETSC_COPY_VALUES,&is1);CHKERRQ(ierr);
498a7e14dcfSSatish Balay     ierr = VecScatterCreate(ipmP->bigstep,sis,ipmP->s,is1,&ipmP->step2);CHKERRQ(ierr);
499a7e14dcfSSatish Balay     ierr = ISDestroy(&sis);CHKERRQ(ierr);
500a7e14dcfSSatish Balay 
501a7e14dcfSSatish Balay     for (i=sstart;i<send;i++) {
502a7e14dcfSSatish Balay       stepind[i-sstart] = i+ipmP->n+ipmP->me;
503a7e14dcfSSatish Balay       cind[i-sstart] = i;
504a7e14dcfSSatish Balay     }
505a7e14dcfSSatish Balay     ierr = ISCreateGeneral(comm,send-sstart,stepind,PETSC_COPY_VALUES,&sis);CHKERRQ(ierr);
506a7e14dcfSSatish Balay     ierr = VecScatterCreate(ipmP->s,is1,ipmP->bigrhs,sis,&ipmP->rhs3);CHKERRQ(ierr);
507a7e14dcfSSatish Balay     ierr = ISDestroy(&sis);CHKERRQ(ierr);
508a7e14dcfSSatish Balay     ierr = ISDestroy(&is1);CHKERRQ(ierr);
509a7e14dcfSSatish Balay   }
510a7e14dcfSSatish Balay 
511a7e14dcfSSatish Balay   if (ipmP->me > 0) {
512a7e14dcfSSatish Balay     ierr = VecGetOwnershipRange(tao->constraints_equality,&ucestart,&uceend);CHKERRQ(ierr);
513a7e14dcfSSatish Balay     for (i=ucestart;i<uceend;i++) {
514a7e14dcfSSatish Balay       stepind[i-ucestart] = i + ipmP->n+ipmP->nb;
515a7e14dcfSSatish Balay       uceind[i-ucestart] = i;
516a7e14dcfSSatish Balay     }
517a7e14dcfSSatish Balay 
518a7e14dcfSSatish Balay     ierr = ISCreateGeneral(comm,uceend-ucestart,stepind,PETSC_COPY_VALUES,&sis);CHKERRQ(ierr);
519a7e14dcfSSatish Balay     ierr = ISCreateGeneral(comm,uceend-ucestart,uceind,PETSC_COPY_VALUES,&is1);CHKERRQ(ierr);
520a7e14dcfSSatish Balay     ierr = VecScatterCreate(ipmP->bigstep,sis,tao->constraints_equality,is1,&ipmP->step3);CHKERRQ(ierr);
521a7e14dcfSSatish Balay     ierr = ISDestroy(&sis);CHKERRQ(ierr);
522a7e14dcfSSatish Balay 
523a7e14dcfSSatish Balay 
524a7e14dcfSSatish Balay     for (i=ucestart;i<uceend;i++) {
525a7e14dcfSSatish Balay       stepind[i-ucestart] = i + ipmP->n;
526a7e14dcfSSatish Balay     }
527a7e14dcfSSatish Balay 
528a7e14dcfSSatish Balay     ierr = ISCreateGeneral(comm,uceend-ucestart,stepind,PETSC_COPY_VALUES,&sis);CHKERRQ(ierr);
529a7e14dcfSSatish Balay     ierr = VecScatterCreate(tao->constraints_equality,is1,ipmP->bigrhs,sis,&ipmP->rhs2);CHKERRQ(ierr);
530a7e14dcfSSatish Balay     ierr = ISDestroy(&sis);CHKERRQ(ierr);
531a7e14dcfSSatish Balay     ierr = ISDestroy(&is1);CHKERRQ(ierr);
532a7e14dcfSSatish Balay   }
533a7e14dcfSSatish Balay 
534a7e14dcfSSatish Balay   if (ipmP->nb > 0) {
535a7e14dcfSSatish Balay     for (i=sstart;i<send;i++) {
536a7e14dcfSSatish Balay       stepind[i-sstart] = i + ipmP->n + ipmP->nb + ipmP->me;
537a7e14dcfSSatish Balay       cind[i-sstart] = i;
538a7e14dcfSSatish Balay     }
539a7e14dcfSSatish Balay     ierr = ISCreateGeneral(comm,send-sstart,cind,PETSC_COPY_VALUES,&is1);
540a7e14dcfSSatish Balay     ierr = ISCreateGeneral(comm,send-sstart,stepind,PETSC_COPY_VALUES,&sis);CHKERRQ(ierr);
541a7e14dcfSSatish Balay     ierr = VecScatterCreate(ipmP->bigstep,sis,ipmP->s,is1,&ipmP->step4);CHKERRQ(ierr);
542a7e14dcfSSatish Balay     ierr = VecScatterCreate(ipmP->s,is1,ipmP->bigrhs,sis,&ipmP->rhs4);CHKERRQ(ierr);
543a7e14dcfSSatish Balay     ierr = ISDestroy(&sis);CHKERRQ(ierr);
544a7e14dcfSSatish Balay     ierr = ISDestroy(&is1);CHKERRQ(ierr);
545a7e14dcfSSatish Balay   }
546a7e14dcfSSatish Balay 
547a7e14dcfSSatish Balay   ierr = PetscFree(stepind);CHKERRQ(ierr);
548a7e14dcfSSatish Balay   ierr = PetscFree(cind);CHKERRQ(ierr);
549a7e14dcfSSatish Balay   ierr = PetscFree(ucind);CHKERRQ(ierr);
550a7e14dcfSSatish Balay   ierr = PetscFree(uceind);CHKERRQ(ierr);
551a7e14dcfSSatish Balay   ierr = PetscFree(xind);CHKERRQ(ierr);
552a7e14dcfSSatish Balay   PetscFunctionReturn(0);
553a7e14dcfSSatish Balay }
554a7e14dcfSSatish Balay 
555a7e14dcfSSatish Balay #undef __FUNCT__
556a7e14dcfSSatish Balay #define __FUNCT__ "TaoDestroy_IPM"
557a7e14dcfSSatish Balay static PetscErrorCode TaoDestroy_IPM(TaoSolver tao)
558a7e14dcfSSatish Balay {
559a7e14dcfSSatish Balay   TAO_IPM        *ipmP = (TAO_IPM*)tao->data;
560a7e14dcfSSatish Balay   PetscErrorCode ierr;
56147a47007SBarry Smith 
562a7e14dcfSSatish Balay   PetscFunctionBegin;
563a7e14dcfSSatish Balay   ierr = VecDestroy(&ipmP->rd);CHKERRQ(ierr);
564a7e14dcfSSatish Balay   ierr = VecDestroy(&ipmP->rpe);CHKERRQ(ierr);
565a7e14dcfSSatish Balay   ierr = VecDestroy(&ipmP->rpi);CHKERRQ(ierr);
566a7e14dcfSSatish Balay   ierr = VecDestroy(&ipmP->work);CHKERRQ(ierr);
567a7e14dcfSSatish Balay   ierr = VecDestroy(&ipmP->lamdae);CHKERRQ(ierr);
568a7e14dcfSSatish Balay   ierr = VecDestroy(&ipmP->lamdai);CHKERRQ(ierr);
569a7e14dcfSSatish Balay   ierr = VecDestroy(&ipmP->s);CHKERRQ(ierr);
570a7e14dcfSSatish Balay   ierr = VecDestroy(&ipmP->ds);CHKERRQ(ierr);
571a7e14dcfSSatish Balay   ierr = VecDestroy(&ipmP->ci);CHKERRQ(ierr);
572a7e14dcfSSatish Balay 
573a7e14dcfSSatish Balay   ierr = VecDestroy(&ipmP->rhs_x);CHKERRQ(ierr);
574a7e14dcfSSatish Balay   ierr = VecDestroy(&ipmP->rhs_lamdae);CHKERRQ(ierr);
575a7e14dcfSSatish Balay   ierr = VecDestroy(&ipmP->rhs_lamdai);CHKERRQ(ierr);
576a7e14dcfSSatish Balay   ierr = VecDestroy(&ipmP->rhs_s);CHKERRQ(ierr);
577a7e14dcfSSatish Balay 
578a7e14dcfSSatish Balay   ierr = VecDestroy(&ipmP->save_x);CHKERRQ(ierr);
579a7e14dcfSSatish Balay   ierr = VecDestroy(&ipmP->save_lamdae);CHKERRQ(ierr);
580a7e14dcfSSatish Balay   ierr = VecDestroy(&ipmP->save_lamdai);CHKERRQ(ierr);
581a7e14dcfSSatish Balay   ierr = VecDestroy(&ipmP->save_s);CHKERRQ(ierr);
582a7e14dcfSSatish Balay 
583a7e14dcfSSatish Balay   ierr = VecScatterDestroy(&ipmP->step1);CHKERRQ(ierr);
584a7e14dcfSSatish Balay   ierr = VecScatterDestroy(&ipmP->step2);CHKERRQ(ierr);
585a7e14dcfSSatish Balay   ierr = VecScatterDestroy(&ipmP->step3);CHKERRQ(ierr);
586a7e14dcfSSatish Balay   ierr = VecScatterDestroy(&ipmP->step4);CHKERRQ(ierr);
587a7e14dcfSSatish Balay 
588a7e14dcfSSatish Balay   ierr = VecScatterDestroy(&ipmP->rhs1);CHKERRQ(ierr);
589a7e14dcfSSatish Balay   ierr = VecScatterDestroy(&ipmP->rhs2);CHKERRQ(ierr);
590a7e14dcfSSatish Balay   ierr = VecScatterDestroy(&ipmP->rhs3);CHKERRQ(ierr);
591a7e14dcfSSatish Balay   ierr = VecScatterDestroy(&ipmP->rhs4);CHKERRQ(ierr);
592a7e14dcfSSatish Balay 
593a7e14dcfSSatish Balay   ierr = VecScatterDestroy(&ipmP->ci_scat);CHKERRQ(ierr);
594a7e14dcfSSatish Balay   ierr = VecScatterDestroy(&ipmP->xl_scat);CHKERRQ(ierr);
595a7e14dcfSSatish Balay   ierr = VecScatterDestroy(&ipmP->xu_scat);CHKERRQ(ierr);
596a7e14dcfSSatish Balay 
597a7e14dcfSSatish Balay   ierr = VecDestroy(&ipmP->dlamdai);CHKERRQ(ierr);
598a7e14dcfSSatish Balay   ierr = VecDestroy(&ipmP->dlamdae);CHKERRQ(ierr);
599a7e14dcfSSatish Balay   ierr = VecDestroy(&ipmP->Zero_nb);CHKERRQ(ierr);
600a7e14dcfSSatish Balay   ierr = VecDestroy(&ipmP->One_nb);CHKERRQ(ierr);
601a7e14dcfSSatish Balay   ierr = VecDestroy(&ipmP->Inf_nb);CHKERRQ(ierr);
602a7e14dcfSSatish Balay   ierr = VecDestroy(&ipmP->complementarity);CHKERRQ(ierr);
603a7e14dcfSSatish Balay 
604a7e14dcfSSatish Balay   ierr = VecDestroy(&ipmP->bigrhs);CHKERRQ(ierr);
605a7e14dcfSSatish Balay   ierr = VecDestroy(&ipmP->bigstep);CHKERRQ(ierr);
606a7e14dcfSSatish Balay   ierr = MatDestroy(&ipmP->Ai);CHKERRQ(ierr);
607a7e14dcfSSatish Balay   ierr = MatDestroy(&ipmP->K);CHKERRQ(ierr);
608a7e14dcfSSatish Balay   ierr = ISDestroy(&ipmP->isxu);CHKERRQ(ierr);
609a7e14dcfSSatish Balay   ierr = ISDestroy(&ipmP->isxl);CHKERRQ(ierr);
610a7e14dcfSSatish Balay   ierr = PetscFree(tao->data);CHKERRQ(ierr);
6116c23d075SBarry Smith   tao->data = NULL;
612a7e14dcfSSatish Balay   PetscFunctionReturn(0);
613a7e14dcfSSatish Balay }
614a7e14dcfSSatish Balay 
615a7e14dcfSSatish Balay #undef __FUNCT__
616a7e14dcfSSatish Balay #define __FUNCT__ "TaoSetFromOptions_IPM"
617a7e14dcfSSatish Balay static PetscErrorCode TaoSetFromOptions_IPM(TaoSolver tao)
618a7e14dcfSSatish Balay {
619a7e14dcfSSatish Balay   TAO_IPM        *ipmP = (TAO_IPM*)tao->data;
620a7e14dcfSSatish Balay   PetscErrorCode ierr;
621a7e14dcfSSatish Balay   PetscBool      flg;
62247a47007SBarry Smith 
623a7e14dcfSSatish Balay   PetscFunctionBegin;
624a7e14dcfSSatish Balay   ierr = PetscOptionsHead("IPM method for constrained optimization");CHKERRQ(ierr);
6256c23d075SBarry Smith   ierr = PetscOptionsBool("-ipm_monitorkkt","monitor kkt status",NULL,ipmP->monitorkkt,&ipmP->monitorkkt,&flg);CHKERRQ(ierr);
6266c23d075SBarry Smith   ierr = PetscOptionsReal("-ipm_pushs","parameter to push initial slack variables away from bounds",NULL,ipmP->pushs,&ipmP->pushs,&flg);
6276c23d075SBarry Smith   ierr = PetscOptionsReal("-ipm_pushnu","parameter to push initial (inequality) dual variables away from bounds",NULL,ipmP->pushnu,&ipmP->pushnu,&flg);
628a7e14dcfSSatish Balay   ierr = PetscOptionsTail();CHKERRQ(ierr);
629a7e14dcfSSatish Balay   ierr =KSPSetFromOptions(tao->ksp);CHKERRQ(ierr);
630a7e14dcfSSatish Balay   PetscFunctionReturn(0);
631a7e14dcfSSatish Balay }
632a7e14dcfSSatish Balay 
633a7e14dcfSSatish Balay #undef __FUNCT__
634a7e14dcfSSatish Balay #define __FUNCT__ "TaoView_IPM"
635a7e14dcfSSatish Balay static PetscErrorCode TaoView_IPM(TaoSolver tao, PetscViewer viewer)
636a7e14dcfSSatish Balay {
637a7e14dcfSSatish Balay   return 0;
638a7e14dcfSSatish Balay }
639a7e14dcfSSatish Balay 
640a7e14dcfSSatish Balay /* IPMObjectiveAndGradient()
641a7e14dcfSSatish Balay    f = d'x + 0.5 * x' * H * x
642a7e14dcfSSatish Balay    rd = H*x + d + Ae'*lame - Ai'*lami
643a7e14dcfSSatish Balay    rpe = Ae*x - be
644a7e14dcfSSatish Balay    rpi = Ai*x - yi - bi
645a7e14dcfSSatish Balay    mu = yi' * lami/mi;
646a7e14dcfSSatish Balay    com = yi.*lami
647a7e14dcfSSatish Balay 
648a7e14dcfSSatish Balay    phi = ||rd|| + ||rpe|| + ||rpi|| + ||com||
649a7e14dcfSSatish Balay */
650a7e14dcfSSatish Balay /*
651a7e14dcfSSatish Balay #undef __FUNCT__
652a7e14dcfSSatish Balay #define __FUNCT__ "IPMObjective"
653a7e14dcfSSatish Balay static PetscErrorCode IPMObjective(TaoLineSearch ls, Vec X, PetscReal *f, void *tptr)
654a7e14dcfSSatish Balay {
655a7e14dcfSSatish Balay   TaoSolver tao = (TaoSolver)tptr;
656a7e14dcfSSatish Balay   TAO_IPM *ipmP = (TAO_IPM*)tao->data;
657a7e14dcfSSatish Balay   PetscErrorCode ierr;
658a7e14dcfSSatish Balay   PetscFunctionBegin;
659a7e14dcfSSatish Balay   ierr = IPMComputeKKT(tao);CHKERRQ(ierr);
660a7e14dcfSSatish Balay   *f = ipmP->phi;
661a7e14dcfSSatish Balay   PetscFunctionReturn(0);
662a7e14dcfSSatish Balay }
663a7e14dcfSSatish Balay */
664a7e14dcfSSatish Balay 
665a7e14dcfSSatish Balay /*
666a7e14dcfSSatish Balay    f = d'x + 0.5 * x' * H * x
667a7e14dcfSSatish Balay    rd = H*x + d + Ae'*lame - Ai'*lami
668a7e14dcfSSatish Balay        Ai =   jac_ineq
669a7e14dcfSSatish Balay                I (w/lb)
670a7e14dcfSSatish Balay               -I (w/ub)
671a7e14dcfSSatish Balay 
672a7e14dcfSSatish Balay    rpe = ce
673a7e14dcfSSatish Balay    rpi = ci - s;
674a7e14dcfSSatish Balay    com = s.*lami
675a7e14dcfSSatish Balay    mu = yi' * lami/mi;
676a7e14dcfSSatish Balay 
677a7e14dcfSSatish Balay    phi = ||rd|| + ||rpe|| + ||rpi|| + ||com||
678a7e14dcfSSatish Balay */
679a7e14dcfSSatish Balay #undef __FUNCT__
680a7e14dcfSSatish Balay #define __FUNCT__ "IPMComputeKKT"
681a7e14dcfSSatish Balay static PetscErrorCode IPMComputeKKT(TaoSolver tao)
682a7e14dcfSSatish Balay {
683a7e14dcfSSatish Balay   TAO_IPM        *ipmP = (TAO_IPM *)tao->data;
684a7e14dcfSSatish Balay   PetscScalar    norm;
685a7e14dcfSSatish Balay   PetscErrorCode ierr;
68647a47007SBarry Smith 
68747a47007SBarry Smith   PetscFunctionBegin;
688a7e14dcfSSatish Balay   ierr = VecCopy(tao->gradient,ipmP->rd);CHKERRQ(ierr);
689a7e14dcfSSatish Balay 
690a7e14dcfSSatish Balay   if (ipmP->me > 0) {
691a7e14dcfSSatish Balay     /* rd = gradient + Ae'*lamdae */
692a7e14dcfSSatish Balay     ierr = MatMultTranspose(tao->jacobian_equality,ipmP->lamdae,ipmP->work);CHKERRQ(ierr);
693a7e14dcfSSatish Balay     ierr = VecAXPY(ipmP->rd, 1.0, ipmP->work);CHKERRQ(ierr);
694a7e14dcfSSatish Balay 
695a7e14dcfSSatish Balay #if defined DEBUG_KKT
696a7e14dcfSSatish Balay     PetscPrintf(PETSC_COMM_WORLD,"\nAe.lamdae\n");
697a7e14dcfSSatish Balay     ierr = VecView(ipmP->work,0);
698a7e14dcfSSatish Balay #endif
699a7e14dcfSSatish Balay     /* rpe = ce(x) */
700a7e14dcfSSatish Balay     ierr = VecCopy(tao->constraints_equality,ipmP->rpe);CHKERRQ(ierr);
701a7e14dcfSSatish Balay 
702a7e14dcfSSatish Balay   }
703a7e14dcfSSatish Balay   if (ipmP->nb > 0) {
704a7e14dcfSSatish Balay     /* rd = rd - Ai'*lamdai */
705a7e14dcfSSatish Balay     ierr = MatMultTranspose(ipmP->Ai,ipmP->lamdai,ipmP->work);CHKERRQ(ierr);
706a7e14dcfSSatish Balay     ierr = VecAXPY(ipmP->rd, -1.0, ipmP->work);CHKERRQ(ierr);
707a7e14dcfSSatish Balay #if defined DEBUG_KKT
708a7e14dcfSSatish Balay     PetscPrintf(PETSC_COMM_WORLD,"\nAi\n");
709a7e14dcfSSatish Balay     ierr = MatView(ipmP->Ai,0);
710a7e14dcfSSatish Balay     PetscPrintf(PETSC_COMM_WORLD,"\nAi.lamdai\n");
711a7e14dcfSSatish Balay     ierr = VecView(ipmP->work,0);
712a7e14dcfSSatish Balay #endif
713a7e14dcfSSatish Balay     /* rpi = cin - s */
714a7e14dcfSSatish Balay     ierr = VecCopy(ipmP->ci,ipmP->rpi);CHKERRQ(ierr);
715a7e14dcfSSatish Balay     ierr = VecAXPY(ipmP->rpi, -1.0, ipmP->s);CHKERRQ(ierr);
716a7e14dcfSSatish Balay 
717a7e14dcfSSatish Balay     /* com = s .* lami */
718a7e14dcfSSatish Balay     ierr = VecPointwiseMult(ipmP->complementarity, ipmP->s,ipmP->lamdai);CHKERRQ(ierr);
719a7e14dcfSSatish Balay 
720a7e14dcfSSatish Balay   }
721a7e14dcfSSatish Balay   /* phi = ||rd; rpe; rpi; com|| */
722a7e14dcfSSatish Balay   ierr = VecDot(ipmP->rd,ipmP->rd,&norm);CHKERRQ(ierr);
723a7e14dcfSSatish Balay   ipmP->phi = norm;
724a7e14dcfSSatish Balay   if (ipmP->me > 0 ) {
725a7e14dcfSSatish Balay     ierr = VecDot(ipmP->rpe,ipmP->rpe,&norm);CHKERRQ(ierr);
726a7e14dcfSSatish Balay     ipmP->phi += norm;
727a7e14dcfSSatish Balay   }
728a7e14dcfSSatish Balay   if (ipmP->nb > 0) {
729a7e14dcfSSatish Balay     ierr = VecDot(ipmP->rpi,ipmP->rpi,&norm);CHKERRQ(ierr);
730a7e14dcfSSatish Balay     ipmP->phi += norm;
731a7e14dcfSSatish Balay     ierr = VecDot(ipmP->complementarity,ipmP->complementarity,&norm);CHKERRQ(ierr);
732a7e14dcfSSatish Balay     ipmP->phi += norm;
733a7e14dcfSSatish Balay     /* mu = s'*lami/nb */
734a7e14dcfSSatish Balay     ierr = VecDot(ipmP->s,ipmP->lamdai,&ipmP->mu);CHKERRQ(ierr);
735a7e14dcfSSatish Balay     ipmP->mu /= ipmP->nb;
736a7e14dcfSSatish Balay   } else {
737a7e14dcfSSatish Balay     ipmP->mu = 1.0;
738a7e14dcfSSatish Balay   }
739a7e14dcfSSatish Balay 
740a7e14dcfSSatish Balay   ipmP->phi = PetscSqrtScalar(ipmP->phi);
741a7e14dcfSSatish Balay #if defined DEBUG_KKT
742a7e14dcfSSatish Balay   if (ipmP->monitorkkt) {
74347a47007SBarry Smith     ierr = PetscPrintf(PETSC_COMM_WORLD,"obj=%g,\tphi = %g,\tmu=%g\talpha1=%g\talpha2=%g\n",(double)ipmP->kkt_f,(double)ipmP->phi,(double)ipmP->mu,(double)ipmP->alpha1,(double)ipmP->alpha2);
744a7e14dcfSSatish Balay   }
745a7e14dcfSSatish Balay   PetscPrintf(PETSC_COMM_WORLD,"\ngradient\n");
746a7e14dcfSSatish Balay   ierr = VecView(tao->gradient,0);
747a7e14dcfSSatish Balay   PetscPrintf(PETSC_COMM_WORLD,"\nrd\n");
748a7e14dcfSSatish Balay   ierr = VecView(ipmP->rd,0);
749a7e14dcfSSatish Balay   PetscPrintf(PETSC_COMM_WORLD,"\nrpe\n");
750a7e14dcfSSatish Balay   ierr = VecView(ipmP->rpe,0);
751a7e14dcfSSatish Balay   PetscPrintf(PETSC_COMM_WORLD,"\nrpi\n");
752a7e14dcfSSatish Balay   ierr = VecView(ipmP->rpi,0);
753a7e14dcfSSatish Balay   PetscPrintf(PETSC_COMM_WORLD,"\ncomplementarity\n");
754a7e14dcfSSatish Balay   ierr = VecView(ipmP->complementarity,0);
755a7e14dcfSSatish Balay #endif
756a7e14dcfSSatish Balay   PetscFunctionReturn(0);
757a7e14dcfSSatish Balay }
758a7e14dcfSSatish Balay 
759a7e14dcfSSatish Balay #undef __FUNCT__
760a7e14dcfSSatish Balay #define __FUNCT__ "IPMEvaluate"
761a7e14dcfSSatish Balay /* evaluate user info at current point */
762a7e14dcfSSatish Balay PetscErrorCode IPMEvaluate(TaoSolver tao)
763a7e14dcfSSatish Balay {
764a7e14dcfSSatish Balay   TAO_IPM        *ipmP = (TAO_IPM *)tao->data;
765a7e14dcfSSatish Balay   PetscErrorCode ierr;
76647a47007SBarry Smith 
767a7e14dcfSSatish Balay   PetscFunctionBegin;
768a7e14dcfSSatish Balay   ierr = TaoComputeObjectiveAndGradient(tao,tao->solution,&ipmP->kkt_f,tao->gradient);CHKERRQ(ierr);
769a7e14dcfSSatish Balay   ierr = TaoComputeHessian(tao,tao->solution,&tao->hessian,&tao->hessian_pre,&ipmP->Hflag);CHKERRQ(ierr);
770a7e14dcfSSatish Balay   if (ipmP->me > 0) {
771a7e14dcfSSatish Balay     ierr = TaoComputeEqualityConstraints(tao,tao->solution,tao->constraints_equality);
772a7e14dcfSSatish Balay     ierr = TaoComputeJacobianEquality(tao,tao->solution,&tao->jacobian_equality,&tao->jacobian_equality_pre,&ipmP->Aiflag);CHKERRQ(ierr);
773a7e14dcfSSatish Balay   }
774a7e14dcfSSatish Balay   if (ipmP->mi > 0) {
775a7e14dcfSSatish Balay     ierr = TaoComputeInequalityConstraints(tao,tao->solution,tao->constraints_inequality);
776a7e14dcfSSatish Balay     ierr = TaoComputeJacobianInequality(tao,tao->solution,&tao->jacobian_inequality,&tao->jacobian_inequality_pre,&ipmP->Aeflag);CHKERRQ(ierr);
777a7e14dcfSSatish Balay 
778a7e14dcfSSatish Balay   }
779a7e14dcfSSatish Balay   if (ipmP->nb > 0) {
780a7e14dcfSSatish Balay     /* Ai' =   jac_ineq | I (w/lb) | -I (w/ub)  */
781a7e14dcfSSatish Balay     ierr = IPMUpdateAi(tao);CHKERRQ(ierr);
782a7e14dcfSSatish Balay   }
783a7e14dcfSSatish Balay   PetscFunctionReturn(0);
784a7e14dcfSSatish Balay }
785a7e14dcfSSatish Balay 
786a7e14dcfSSatish Balay #undef __FUNCT__
787a7e14dcfSSatish Balay #define __FUNCT__ "IPMPushInitialPoint"
788a7e14dcfSSatish Balay /* Push initial point away from bounds */
789a7e14dcfSSatish Balay PetscErrorCode IPMPushInitialPoint(TaoSolver tao)
790a7e14dcfSSatish Balay {
791a7e14dcfSSatish Balay   TAO_IPM        *ipmP = (TAO_IPM *)tao->data;
792a7e14dcfSSatish Balay   PetscErrorCode ierr;
793a7e14dcfSSatish Balay 
79447a47007SBarry Smith   PetscFunctionBegin;
795a7e14dcfSSatish Balay   ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr);
796a7e14dcfSSatish Balay   if (tao->XL && tao->XU) {
797a7e14dcfSSatish Balay     ierr = VecMedian(tao->XL, tao->solution, tao->XU, tao->solution);CHKERRQ(ierr);
798a7e14dcfSSatish Balay   }
799a7e14dcfSSatish Balay   if (ipmP->nb > 0) {
800a7e14dcfSSatish Balay     ierr = VecSet(ipmP->s,ipmP->pushs);CHKERRQ(ierr);
801a7e14dcfSSatish Balay     ierr = VecSet(ipmP->lamdai,ipmP->pushnu);CHKERRQ(ierr);
802a7e14dcfSSatish Balay     if (ipmP->mi > 0) {
803a7e14dcfSSatish Balay       ierr = VecSet(tao->DI,ipmP->pushnu);CHKERRQ(ierr);
804a7e14dcfSSatish Balay     }
805a7e14dcfSSatish Balay   }
806a7e14dcfSSatish Balay   if (ipmP->me > 0) {
807a7e14dcfSSatish Balay     ierr = VecSet(tao->DE,1.0);CHKERRQ(ierr);
808a7e14dcfSSatish Balay     ierr = VecSet(ipmP->lamdae,1.0);CHKERRQ(ierr);
809a7e14dcfSSatish Balay   }
810a7e14dcfSSatish Balay   PetscFunctionReturn(0);
811a7e14dcfSSatish Balay }
812a7e14dcfSSatish Balay 
813a7e14dcfSSatish Balay #undef __FUNCT__
814a7e14dcfSSatish Balay #define __FUNCT__ "IPMUpdateAi"
815a7e14dcfSSatish Balay PetscErrorCode IPMUpdateAi(TaoSolver tao)
816a7e14dcfSSatish Balay {
817a7e14dcfSSatish Balay   /* Ai =     Ji
818a7e14dcfSSatish Balay               I (w/lb)
819a7e14dcfSSatish Balay              -I (w/ub) */
820a7e14dcfSSatish Balay 
821a7e14dcfSSatish Balay   /* Ci =    user->ci
822a7e14dcfSSatish Balay              Xi - lb (w/lb)
823a7e14dcfSSatish Balay              -Xi + ub (w/ub)  */
824a7e14dcfSSatish Balay 
825a7e14dcfSSatish Balay   TAO_IPM           *ipmP = (TAO_IPM *)tao->data;
826a7e14dcfSSatish Balay   MPI_Comm          comm;
827a7e14dcfSSatish Balay   PetscInt          i;
828a7e14dcfSSatish Balay   PetscScalar       newval;
829a7e14dcfSSatish Balay   PetscInt          newrow,newcol,ncols;
830a7e14dcfSSatish Balay   const PetscScalar *vals;
831a7e14dcfSSatish Balay   const PetscInt    *cols;
832a7e14dcfSSatish Balay   PetscInt          astart,aend,jstart,jend;
833a7e14dcfSSatish Balay   PetscInt          *nonzeros;
834a7e14dcfSSatish Balay   PetscInt          r2,r3,r4;
83547a47007SBarry Smith   PetscMPIInt       size;
836a7e14dcfSSatish Balay   PetscErrorCode    ierr;
837a7e14dcfSSatish Balay 
838a7e14dcfSSatish Balay   PetscFunctionBegin;
839a7e14dcfSSatish Balay   CHKMEMQ;
840a7e14dcfSSatish Balay   r2 = ipmP->mi;
841a7e14dcfSSatish Balay   r3 = r2 + ipmP->nxlb;
842a7e14dcfSSatish Balay   r4 = r3 + ipmP->nxub;
843a7e14dcfSSatish Balay 
844a7e14dcfSSatish Balay   if (!ipmP->nb) {
845a7e14dcfSSatish Balay     PetscFunctionReturn(0);
846a7e14dcfSSatish Balay   }
847a7e14dcfSSatish Balay   CHKMEMQ;
848a7e14dcfSSatish Balay 
849a7e14dcfSSatish Balay   /* Create Ai matrix if it doesn't exist yet */
850a7e14dcfSSatish Balay   if (!ipmP->Ai) {
851a7e14dcfSSatish Balay     comm = ((PetscObject)(tao->solution))->comm;
8520e660641SBarry Smith     ierr = PetscMalloc1(ipmP->nb,&nonzeros);CHKERRQ(ierr);
85347a47007SBarry Smith     ierr = MPI_Comm_size(comm,&size);
85447a47007SBarry Smith     if (size == 1) {
855a7e14dcfSSatish Balay       for (i=0;i<ipmP->mi;i++) {
8566c23d075SBarry Smith         ierr = MatGetRow(tao->jacobian_inequality,i,&ncols,NULL,NULL);CHKERRQ(ierr);
857a7e14dcfSSatish Balay         nonzeros[i] = ncols;
8586c23d075SBarry Smith         ierr = MatRestoreRow(tao->jacobian_inequality,i,&ncols,NULL,NULL);CHKERRQ(ierr);
859a7e14dcfSSatish Balay       }
860a7e14dcfSSatish Balay       for (i=r2;i<r4;i++) {
861a7e14dcfSSatish Balay         nonzeros[i] = 1;
862a7e14dcfSSatish Balay       }
863a7e14dcfSSatish Balay     }
864a7e14dcfSSatish Balay     ierr = MatCreate(comm,&ipmP->Ai);CHKERRQ(ierr);
865a7e14dcfSSatish Balay     ierr = MatSetType(ipmP->Ai,MATAIJ);CHKERRQ(ierr);
866a7e14dcfSSatish Balay     ierr = MatSetSizes(ipmP->Ai,PETSC_DECIDE,PETSC_DECIDE,ipmP->nb,ipmP->n);CHKERRQ(ierr);
867a7e14dcfSSatish Balay     ierr = MatSetFromOptions(ipmP->Ai);CHKERRQ(ierr);
8686c23d075SBarry Smith     ierr = MatMPIAIJSetPreallocation(ipmP->Ai,ipmP->nb,NULL,ipmP->nb,NULL);
869a7e14dcfSSatish Balay     ierr = MatSeqAIJSetPreallocation(ipmP->Ai,PETSC_DEFAULT,nonzeros);CHKERRQ(ierr);
87047a47007SBarry Smith     if (size ==1) {
871a7e14dcfSSatish Balay       ierr = PetscFree(nonzeros);CHKERRQ(ierr);
872a7e14dcfSSatish Balay     }
873a7e14dcfSSatish Balay   }
874a7e14dcfSSatish Balay 
875a7e14dcfSSatish Balay 
876a7e14dcfSSatish Balay   /* Copy values from user jacobian to Ai */
877a7e14dcfSSatish Balay   ierr = MatGetOwnershipRange(ipmP->Ai,&astart,&aend);CHKERRQ(ierr);
878a7e14dcfSSatish Balay 
879a7e14dcfSSatish Balay   /* Ai w/lb */
880a7e14dcfSSatish Balay   if (ipmP->mi) {
881a7e14dcfSSatish Balay     ierr = MatZeroEntries(ipmP->Ai);CHKERRQ(ierr);
882a7e14dcfSSatish Balay     ierr = MatGetOwnershipRange(tao->jacobian_inequality,&jstart,&jend);CHKERRQ(ierr);
883a7e14dcfSSatish Balay     for (i=jstart;i<jend;i++) {
884a7e14dcfSSatish Balay       ierr = MatGetRow(tao->jacobian_inequality,i,&ncols,&cols,&vals);CHKERRQ(ierr);
885a7e14dcfSSatish Balay       newrow = i;
886a7e14dcfSSatish Balay       ierr = MatSetValues(ipmP->Ai,1,&newrow,ncols,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
887a7e14dcfSSatish Balay       ierr = MatRestoreRow(tao->jacobian_inequality,i,&ncols,&cols,&vals);CHKERRQ(ierr);
888a7e14dcfSSatish Balay     }
889a7e14dcfSSatish Balay   }
890a7e14dcfSSatish Balay 
891a7e14dcfSSatish Balay 
892a7e14dcfSSatish Balay   /* I w/ xlb */
893a7e14dcfSSatish Balay   if (ipmP->nxlb) {
894a7e14dcfSSatish Balay     for (i=0;i<ipmP->nxlb;i++) {
895a7e14dcfSSatish Balay       if (i>=astart && i<aend) {
896a7e14dcfSSatish Balay         newrow = i+r2;
897a7e14dcfSSatish Balay         newcol = i;
898a7e14dcfSSatish Balay         newval = 1.0;
899a7e14dcfSSatish Balay         ierr = MatSetValues(ipmP->Ai,1,&newrow,1,&newcol,&newval,INSERT_VALUES);CHKERRQ(ierr);
900a7e14dcfSSatish Balay       }
901a7e14dcfSSatish Balay     }
902a7e14dcfSSatish Balay   }
903a7e14dcfSSatish Balay   if (ipmP->nxub) {
904a7e14dcfSSatish Balay     /* I w/ xub */
905a7e14dcfSSatish Balay     for (i=0;i<ipmP->nxub;i++) {
906a7e14dcfSSatish Balay       if (i>=astart && i<aend) {
907a7e14dcfSSatish Balay       newrow = i+r3;
908a7e14dcfSSatish Balay       newcol = i;
909a7e14dcfSSatish Balay       newval = -1.0;
910a7e14dcfSSatish Balay       ierr = MatSetValues(ipmP->Ai,1,&newrow,1,&newcol,&newval,INSERT_VALUES);CHKERRQ(ierr);
911a7e14dcfSSatish Balay       }
912a7e14dcfSSatish Balay     }
913a7e14dcfSSatish Balay   }
914a7e14dcfSSatish Balay 
915a7e14dcfSSatish Balay 
916a7e14dcfSSatish Balay   ierr = MatAssemblyBegin(ipmP->Ai,MAT_FINAL_ASSEMBLY);
917a7e14dcfSSatish Balay   ierr = MatAssemblyEnd(ipmP->Ai,MAT_FINAL_ASSEMBLY);
918a7e14dcfSSatish Balay   CHKMEMQ;
919a7e14dcfSSatish Balay 
920a7e14dcfSSatish Balay   ierr = VecSet(ipmP->ci,0.0);CHKERRQ(ierr);
921a7e14dcfSSatish Balay 
922a7e14dcfSSatish Balay   /* user ci */
923a7e14dcfSSatish Balay   if (ipmP->mi > 0) {
924a7e14dcfSSatish Balay     ierr = VecScatterBegin(ipmP->ci_scat,tao->constraints_inequality,ipmP->ci,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
925a7e14dcfSSatish Balay     ierr = VecScatterEnd(ipmP->ci_scat,tao->constraints_inequality,ipmP->ci,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
926a7e14dcfSSatish Balay   }
927a7e14dcfSSatish Balay   if (!ipmP->work){
928a7e14dcfSSatish Balay     VecDuplicate(tao->solution,&ipmP->work);
929a7e14dcfSSatish Balay   }
930a7e14dcfSSatish Balay   ierr = VecCopy(tao->solution,ipmP->work);CHKERRQ(ierr);
931a7e14dcfSSatish Balay   if (tao->XL) {
932a7e14dcfSSatish Balay     ierr = VecAXPY(ipmP->work,-1.0,tao->XL);CHKERRQ(ierr);
933a7e14dcfSSatish Balay 
934a7e14dcfSSatish Balay     /* lower bounds on variables */
935a7e14dcfSSatish Balay     if (ipmP->nxlb > 0) {
936a7e14dcfSSatish Balay       ierr = VecScatterBegin(ipmP->xl_scat,ipmP->work,ipmP->ci,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
937a7e14dcfSSatish Balay       ierr = VecScatterEnd(ipmP->xl_scat,ipmP->work,ipmP->ci,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
938a7e14dcfSSatish Balay     }
939a7e14dcfSSatish Balay   }
940a7e14dcfSSatish Balay   if (tao->XU) {
941a7e14dcfSSatish Balay     /* upper bounds on variables */
942a7e14dcfSSatish Balay     ierr = VecCopy(tao->solution,ipmP->work);CHKERRQ(ierr);
943a7e14dcfSSatish Balay     ierr = VecScale(ipmP->work,-1.0);CHKERRQ(ierr);
944a7e14dcfSSatish Balay     ierr = VecAXPY(ipmP->work,1.0,tao->XU);CHKERRQ(ierr);
945a7e14dcfSSatish Balay     if (ipmP->nxub > 0) {
946a7e14dcfSSatish Balay       ierr = VecScatterBegin(ipmP->xu_scat,ipmP->work,ipmP->ci,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
947a7e14dcfSSatish Balay       ierr = VecScatterEnd(ipmP->xu_scat,ipmP->work,ipmP->ci,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
948a7e14dcfSSatish Balay     }
949a7e14dcfSSatish Balay   }
950a7e14dcfSSatish Balay 
951a7e14dcfSSatish Balay   PetscFunctionReturn(0);
952a7e14dcfSSatish Balay }
953a7e14dcfSSatish Balay 
954a7e14dcfSSatish Balay #undef __FUNCT__
955a7e14dcfSSatish Balay #define __FUNCT__ "IPMUpdateK"
956a7e14dcfSSatish Balay /* create K = [ Hlag , 0 , Ae', -Ai'];
957a7e14dcfSSatish Balay               [Ae , 0,   0  , 0];
958a7e14dcfSSatish Balay               [Ai ,-I,   0 ,  0];
959a7e14dcfSSatish Balay               [ 0 , S ,  0,   Y ];  */
960a7e14dcfSSatish Balay PetscErrorCode IPMUpdateK(TaoSolver tao)
961a7e14dcfSSatish Balay {
962a7e14dcfSSatish Balay   TAO_IPM         *ipmP = (TAO_IPM *)tao->data;
963a7e14dcfSSatish Balay   MPI_Comm        comm;
96447a47007SBarry Smith   PetscMPIInt     size;
965a7e14dcfSSatish Balay   PetscErrorCode  ierr;
966a7e14dcfSSatish Balay   PetscInt        i,j,row;
967a7e14dcfSSatish Balay   PetscInt        ncols,newcol,newcols[2],newrow;
968a7e14dcfSSatish Balay   const PetscInt  *cols;
969a7e14dcfSSatish Balay   const PetscReal *vals;
970a7e14dcfSSatish Balay   PetscReal       *l,*y;
971a7e14dcfSSatish Balay   PetscReal       *newvals;
972a7e14dcfSSatish Balay   PetscReal       newval;
973a7e14dcfSSatish Balay   PetscInt        subsize;
974a7e14dcfSSatish Balay   const PetscInt  *indices;
975a7e14dcfSSatish Balay   PetscInt        *nonzeros,*d_nonzeros,*o_nonzeros;
976a7e14dcfSSatish Balay   PetscInt        bigsize;
977a7e14dcfSSatish Balay   PetscInt        r1,r2,r3;
978a7e14dcfSSatish Balay   PetscInt        c1,c2,c3;
979a7e14dcfSSatish Balay   PetscInt        klocalsize;
980a7e14dcfSSatish Balay   PetscInt        hstart,hend,kstart,kend;
981a7e14dcfSSatish Balay   PetscInt        aistart,aiend,aestart,aeend;
982a7e14dcfSSatish Balay   PetscInt        sstart,send;
983a7e14dcfSSatish Balay 
98447a47007SBarry Smith   PetscFunctionBegin;
985a7e14dcfSSatish Balay   comm = ((PetscObject)(tao->solution))->comm;
98647a47007SBarry Smith   ierr = MPI_Comm_size(comm,&size);
987a7e14dcfSSatish Balay   ierr = IPMUpdateAi(tao);CHKERRQ(ierr);
988a7e14dcfSSatish Balay #if defined DEBUG_K
989a7e14dcfSSatish Balay   PetscPrintf(PETSC_COMM_WORLD,"H\n");  MatView(tao->hessian,0);
990a7e14dcfSSatish Balay   if (ipmP->nb) {
991a7e14dcfSSatish Balay     PetscPrintf(PETSC_COMM_WORLD,"Ai\n"); MatView(ipmP->Ai,0);
992a7e14dcfSSatish Balay   }
993a7e14dcfSSatish Balay   if (ipmP->me) {
994a7e14dcfSSatish Balay     PetscPrintf(PETSC_COMM_WORLD,"Ae\n"); MatView(tao->jacobian_equality,0);
995a7e14dcfSSatish Balay   }
996a7e14dcfSSatish Balay 
997a7e14dcfSSatish Balay #endif
998a7e14dcfSSatish Balay   /* allocate workspace */
999a7e14dcfSSatish Balay   subsize = PetscMax(ipmP->n,ipmP->nb);
1000a7e14dcfSSatish Balay   subsize = PetscMax(ipmP->me,subsize);
1001a7e14dcfSSatish Balay   subsize = PetscMax(2,subsize);
10020e660641SBarry Smith   ierr = PetscMalloc1(subsize,&indices);CHKERRQ(ierr);
10030e660641SBarry Smith   ierr = PetscMalloc1(subsize,&newvals);CHKERRQ(ierr);
1004a7e14dcfSSatish Balay 
1005a7e14dcfSSatish Balay   r1 = c1 = ipmP->n;
1006a7e14dcfSSatish Balay   r2 = r1 + ipmP->me;  c2 = c1 + ipmP->nb;
1007a7e14dcfSSatish Balay   r3 = c3 = r2 + ipmP->nb;
1008a7e14dcfSSatish Balay 
1009a7e14dcfSSatish Balay   bigsize = ipmP->n+2*ipmP->nb+ipmP->me;
1010a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(ipmP->bigrhs,&kstart,&kend);CHKERRQ(ierr);
1011a7e14dcfSSatish Balay   ierr = MatGetOwnershipRange(tao->hessian,&hstart,&hend);CHKERRQ(ierr);
1012a7e14dcfSSatish Balay   klocalsize = kend-kstart;
1013a7e14dcfSSatish Balay   if (!ipmP->K) {
101447a47007SBarry Smith     if (size == 1) {
10150e660641SBarry Smith       ierr = PetscMalloc1((kend-kstart),&nonzeros);CHKERRQ(ierr);
1016a7e14dcfSSatish Balay       for (i=0;i<bigsize;i++) {
1017a7e14dcfSSatish Balay         if (i<r1) {
10186c23d075SBarry Smith           ierr = MatGetRow(tao->hessian,i,&ncols,NULL,NULL);CHKERRQ(ierr);
1019a7e14dcfSSatish Balay           nonzeros[i] = ncols;
10206c23d075SBarry Smith           ierr = MatRestoreRow(tao->hessian,i,&ncols,NULL,NULL);CHKERRQ(ierr);
1021a7e14dcfSSatish Balay           nonzeros[i] += ipmP->me+ipmP->nb;
1022a7e14dcfSSatish Balay         } else if (i<r2) {
1023a7e14dcfSSatish Balay           nonzeros[i-kstart] = ipmP->n;
1024a7e14dcfSSatish Balay         } else if (i<r3) {
1025a7e14dcfSSatish Balay           nonzeros[i-kstart] = ipmP->n+1;
1026a7e14dcfSSatish Balay         } else if (i<bigsize) {
1027a7e14dcfSSatish Balay           nonzeros[i-kstart] = 2;
1028a7e14dcfSSatish Balay         }
1029a7e14dcfSSatish Balay       }
1030a7e14dcfSSatish Balay       ierr = MatCreate(comm,&ipmP->K);CHKERRQ(ierr);
1031a7e14dcfSSatish Balay       ierr = MatSetType(ipmP->K,MATSEQAIJ);CHKERRQ(ierr);
1032a7e14dcfSSatish Balay       ierr = MatSetSizes(ipmP->K,klocalsize,klocalsize,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1033a7e14dcfSSatish Balay       ierr = MatSeqAIJSetPreallocation(ipmP->K,0,nonzeros);CHKERRQ(ierr);
1034a7e14dcfSSatish Balay       ierr = MatSetFromOptions(ipmP->K);CHKERRQ(ierr);
1035a7e14dcfSSatish Balay       ierr = PetscFree(nonzeros);CHKERRQ(ierr);
1036a7e14dcfSSatish Balay     } else {
10370e660641SBarry Smith       ierr = PetscMalloc1((kend-kstart),&d_nonzeros);CHKERRQ(ierr);
10380e660641SBarry Smith       ierr = PetscMalloc1((kend-kstart),&o_nonzeros);CHKERRQ(ierr);
1039a7e14dcfSSatish Balay       for (i=kstart;i<kend;i++) {
1040a7e14dcfSSatish Balay         if (i<r1) {
1041a7e14dcfSSatish Balay           /* TODO fix preallocation for mpi mats */
1042a7e14dcfSSatish Balay           d_nonzeros[i-kstart] = PetscMin(ipmP->n+ipmP->me+ipmP->nb,kend-kstart);
1043a7e14dcfSSatish Balay           o_nonzeros[i-kstart] = PetscMin(ipmP->n+ipmP->me+ipmP->nb,bigsize-(kend-kstart));
1044a7e14dcfSSatish Balay         } else if (i<r2) {
1045a7e14dcfSSatish Balay           d_nonzeros[i-kstart] = PetscMin(ipmP->n,kend-kstart);
1046a7e14dcfSSatish Balay           o_nonzeros[i-kstart] = PetscMin(ipmP->n,bigsize-(kend-kstart));
1047a7e14dcfSSatish Balay         } else if (i<r3) {
1048a7e14dcfSSatish Balay           d_nonzeros[i-kstart] = PetscMin(ipmP->n+2,kend-kstart);
1049a7e14dcfSSatish Balay           o_nonzeros[i-kstart] = PetscMin(ipmP->n+2,bigsize-(kend-kstart));
1050a7e14dcfSSatish Balay         } else {
1051a7e14dcfSSatish Balay           d_nonzeros[i-kstart] = PetscMin(2,kend-kstart);
1052a7e14dcfSSatish Balay           o_nonzeros[i-kstart] = PetscMin(2,bigsize-(kend-kstart));
1053a7e14dcfSSatish Balay         }
1054a7e14dcfSSatish Balay       }
1055a7e14dcfSSatish Balay       ierr = MatCreate(comm,&ipmP->K);CHKERRQ(ierr);
1056a7e14dcfSSatish Balay       ierr = MatSetType(ipmP->K,MATMPIAIJ);CHKERRQ(ierr);
1057a7e14dcfSSatish Balay       ierr = MatSetSizes(ipmP->K,klocalsize,klocalsize,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1058a7e14dcfSSatish Balay       ierr = MatMPIAIJSetPreallocation(ipmP->K,0,d_nonzeros,0,o_nonzeros);CHKERRQ(ierr);
1059a7e14dcfSSatish Balay       ierr = PetscFree(d_nonzeros);CHKERRQ(ierr);
1060a7e14dcfSSatish Balay       ierr = PetscFree(o_nonzeros);CHKERRQ(ierr);
1061a7e14dcfSSatish Balay       ierr = MatSetFromOptions(ipmP->K);CHKERRQ(ierr);
1062a7e14dcfSSatish Balay 
1063a7e14dcfSSatish Balay     }
1064a7e14dcfSSatish Balay   }
1065a7e14dcfSSatish Balay 
1066a7e14dcfSSatish Balay 
1067a7e14dcfSSatish Balay   ierr = MatZeroEntries(ipmP->K);CHKERRQ(ierr);
1068a7e14dcfSSatish Balay   /* Copy H */
1069a7e14dcfSSatish Balay   for (i=hstart;i<hend;i++) {
1070a7e14dcfSSatish Balay     ierr = MatGetRow(tao->hessian,i,&ncols,&cols,&vals);CHKERRQ(ierr);
1071a7e14dcfSSatish Balay     if (ncols > 0) {
1072a7e14dcfSSatish Balay       ierr = MatSetValues(ipmP->K,1,&i,ncols,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
1073a7e14dcfSSatish Balay     }
1074a7e14dcfSSatish Balay     ierr = MatRestoreRow(tao->hessian,i,&ncols,&cols,&vals);CHKERRQ(ierr);
1075a7e14dcfSSatish Balay   }
1076a7e14dcfSSatish Balay 
1077a7e14dcfSSatish Balay   /* Copy Ae and Ae' */
1078a7e14dcfSSatish Balay   if (ipmP->me > 0) {
1079a7e14dcfSSatish Balay     ierr = MatGetOwnershipRange(tao->jacobian_equality,&aestart,&aeend);CHKERRQ(ierr);
1080a7e14dcfSSatish Balay     for (i=aestart;i<aeend;i++) {
1081a7e14dcfSSatish Balay       ierr = MatGetRow(tao->jacobian_equality,i,&ncols,&cols,&vals);CHKERRQ(ierr);
1082a7e14dcfSSatish Balay       if (ncols > 0) {
1083a7e14dcfSSatish Balay         /*Ae*/
1084a7e14dcfSSatish Balay         row = i+r1;
1085a7e14dcfSSatish Balay         ierr = MatSetValues(ipmP->K,1,&row,ncols,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
1086a7e14dcfSSatish Balay         /*Ae'*/
1087a7e14dcfSSatish Balay         for (j=0;j<ncols;j++) {
1088a7e14dcfSSatish Balay           newcol = i + c2;
1089a7e14dcfSSatish Balay           newrow = cols[j];
1090a7e14dcfSSatish Balay           newval = vals[j];
1091a7e14dcfSSatish Balay           ierr = MatSetValues(ipmP->K,1,&newrow,1,&newcol,&newval,INSERT_VALUES);CHKERRQ(ierr);
1092a7e14dcfSSatish Balay         }
1093a7e14dcfSSatish Balay       }
1094a7e14dcfSSatish Balay       ierr = MatRestoreRow(tao->jacobian_equality,i,&ncols,&cols,&vals);CHKERRQ(ierr);
1095a7e14dcfSSatish Balay     }
1096a7e14dcfSSatish Balay   }
1097a7e14dcfSSatish Balay 
1098a7e14dcfSSatish Balay   if (ipmP->nb > 0) {
1099a7e14dcfSSatish Balay     ierr = MatGetOwnershipRange(ipmP->Ai,&aistart,&aiend);CHKERRQ(ierr);
1100a7e14dcfSSatish Balay     /* Copy Ai,and Ai' */
1101a7e14dcfSSatish Balay     for (i=aistart;i<aiend;i++) {
1102a7e14dcfSSatish Balay       row = i+r2;
1103a7e14dcfSSatish Balay       ierr = MatGetRow(ipmP->Ai,i,&ncols,&cols,&vals);CHKERRQ(ierr);
1104a7e14dcfSSatish Balay       if (ncols > 0) {
1105a7e14dcfSSatish Balay         /*Ai*/
1106a7e14dcfSSatish Balay         ierr = MatSetValues(ipmP->K,1,&row,ncols,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
1107a7e14dcfSSatish Balay         /*-Ai'*/
1108a7e14dcfSSatish Balay         for (j=0;j<ncols;j++) {
1109a7e14dcfSSatish Balay           newcol = i + c3;
1110a7e14dcfSSatish Balay           newrow = cols[j];
1111a7e14dcfSSatish Balay           newval = -vals[j];
1112a7e14dcfSSatish Balay           ierr = MatSetValues(ipmP->K,1,&newrow,1,&newcol,&newval,INSERT_VALUES);CHKERRQ(ierr);
1113a7e14dcfSSatish Balay         }
1114a7e14dcfSSatish Balay       }
1115a7e14dcfSSatish Balay       ierr = MatRestoreRow(ipmP->Ai,i,&ncols,&cols,&vals);CHKERRQ(ierr);
1116a7e14dcfSSatish Balay     }
1117a7e14dcfSSatish Balay 
1118a7e14dcfSSatish Balay 
1119a7e14dcfSSatish Balay 
1120a7e14dcfSSatish Balay     /* -I */
1121a7e14dcfSSatish Balay     for (i=kstart;i<kend;i++) {
1122a7e14dcfSSatish Balay       if (i>=r2 && i<r3) {
1123a7e14dcfSSatish Balay         newrow = i;
1124a7e14dcfSSatish Balay         newcol = i-r2+c1;
1125a7e14dcfSSatish Balay         newval = -1.0;
1126a7e14dcfSSatish Balay         MatSetValues(ipmP->K,1,&newrow,1,&newcol,&newval,INSERT_VALUES);CHKERRQ(ierr);
1127a7e14dcfSSatish Balay       }
1128a7e14dcfSSatish Balay     }
1129a7e14dcfSSatish Balay 
1130a7e14dcfSSatish Balay     /* Copy L,Y */
1131a7e14dcfSSatish Balay     ierr = VecGetOwnershipRange(ipmP->s,&sstart,&send);CHKERRQ(ierr);
1132a7e14dcfSSatish Balay     ierr = VecGetArray(ipmP->lamdai,&l);CHKERRQ(ierr);
1133a7e14dcfSSatish Balay     ierr = VecGetArray(ipmP->s,&y);CHKERRQ(ierr);
1134a7e14dcfSSatish Balay 
1135a7e14dcfSSatish Balay     for (i=sstart;i<send;i++) {
1136a7e14dcfSSatish Balay       newcols[0] = c1+i;
1137a7e14dcfSSatish Balay       newcols[1] = c3+i;
1138a7e14dcfSSatish Balay       newvals[0] = l[i-sstart];
1139a7e14dcfSSatish Balay       newvals[1] = y[i-sstart];
1140a7e14dcfSSatish Balay       newrow = r3+i;
1141a7e14dcfSSatish Balay       ierr = MatSetValues(ipmP->K,1,&newrow,2,newcols,newvals,INSERT_VALUES);CHKERRQ(ierr);
1142a7e14dcfSSatish Balay     }
1143a7e14dcfSSatish Balay 
1144a7e14dcfSSatish Balay     ierr = VecRestoreArray(ipmP->lamdai,&l);CHKERRQ(ierr);
1145a7e14dcfSSatish Balay     ierr = VecRestoreArray(ipmP->s,&y);CHKERRQ(ierr);
1146a7e14dcfSSatish Balay   }
1147a7e14dcfSSatish Balay 
1148a7e14dcfSSatish Balay   ierr = PetscFree(indices);CHKERRQ(ierr);
1149a7e14dcfSSatish Balay   ierr = PetscFree(newvals);CHKERRQ(ierr);
1150a7e14dcfSSatish Balay   ierr = MatAssemblyBegin(ipmP->K,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1151a7e14dcfSSatish Balay   ierr = MatAssemblyEnd(ipmP->K,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1152a7e14dcfSSatish Balay #if defined DEBUG_K
1153a7e14dcfSSatish Balay   PetscPrintf(PETSC_COMM_WORLD,"K\n");  MatView(ipmP->K,0);
1154a7e14dcfSSatish Balay #endif
1155a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1156a7e14dcfSSatish Balay }
1157a7e14dcfSSatish Balay 
1158a7e14dcfSSatish Balay #undef __FUNCT__
1159a7e14dcfSSatish Balay #define __FUNCT__ "IPMGatherRHS"
1160a7e14dcfSSatish Balay PetscErrorCode IPMGatherRHS(TaoSolver tao,Vec RHS,Vec X1,Vec X2,Vec X3,Vec X4)
1161a7e14dcfSSatish Balay {
1162a7e14dcfSSatish Balay   TAO_IPM        *ipmP = (TAO_IPM *)tao->data;
1163a7e14dcfSSatish Balay   PetscErrorCode ierr;
1164a7e14dcfSSatish Balay 
116547a47007SBarry Smith   PetscFunctionBegin;
1166a7e14dcfSSatish Balay   /* rhs = [x1      (n)
1167a7e14dcfSSatish Balay             x2     (me)
1168a7e14dcfSSatish Balay             x3     (nb)
1169a7e14dcfSSatish Balay             x4     (nb)] */
1170a7e14dcfSSatish Balay   if (X1) {
1171a7e14dcfSSatish Balay     ierr = VecScatterBegin(ipmP->rhs1,X1,RHS,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1172a7e14dcfSSatish Balay     ierr = VecScatterEnd(ipmP->rhs1,X1,RHS,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1173a7e14dcfSSatish Balay   }
1174a7e14dcfSSatish Balay   if (ipmP->me > 0 && X2) {
1175a7e14dcfSSatish Balay     ierr = VecScatterBegin(ipmP->rhs2,X2,RHS,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1176a7e14dcfSSatish Balay     ierr = VecScatterEnd(ipmP->rhs2,X2,RHS,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1177a7e14dcfSSatish Balay   }
1178a7e14dcfSSatish Balay   if (ipmP->nb > 0) {
1179a7e14dcfSSatish Balay     if (X3) {
1180a7e14dcfSSatish Balay       ierr = VecScatterBegin(ipmP->rhs3,X3,RHS,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1181a7e14dcfSSatish Balay       ierr = VecScatterEnd(ipmP->rhs3,X3,RHS,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1182a7e14dcfSSatish Balay     }
1183a7e14dcfSSatish Balay     if (X4) {
1184a7e14dcfSSatish Balay       ierr = VecScatterBegin(ipmP->rhs4,X4,RHS,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1185a7e14dcfSSatish Balay       ierr = VecScatterEnd(ipmP->rhs4,X4,RHS,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1186a7e14dcfSSatish Balay     }
1187a7e14dcfSSatish Balay   }
1188a7e14dcfSSatish Balay #if defined(DEBUG_SCATTER)
1189a7e14dcfSSatish Balay   PetscPrintf(PETSC_COMM_WORLD,"X1-X4\n");
1190a7e14dcfSSatish Balay   if (X1) {VecView(X1,0);}
1191a7e14dcfSSatish Balay   if (X2) {VecView(X2,0);}
1192a7e14dcfSSatish Balay   if (X3) {VecView(X3,0);}
1193a7e14dcfSSatish Balay   if (X4) {VecView(X4,0);}
1194a7e14dcfSSatish Balay   PetscPrintf(PETSC_COMM_WORLD,"RHS\n");
1195a7e14dcfSSatish Balay   VecView(RHS,0);
1196a7e14dcfSSatish Balay #endif
1197a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1198a7e14dcfSSatish Balay }
1199a7e14dcfSSatish Balay 
1200a7e14dcfSSatish Balay #undef __FUNCT__
1201a7e14dcfSSatish Balay #define __FUNCT__ "IPMScatterStep"
1202a7e14dcfSSatish Balay PetscErrorCode IPMScatterStep(TaoSolver tao, Vec STEP, Vec X1, Vec X2, Vec X3, Vec X4)
1203a7e14dcfSSatish Balay {
1204a7e14dcfSSatish Balay   TAO_IPM        *ipmP = (TAO_IPM *)tao->data;
1205a7e14dcfSSatish Balay   PetscErrorCode ierr;
120647a47007SBarry Smith 
1207a7e14dcfSSatish Balay   PetscFunctionBegin;
1208a7e14dcfSSatish Balay   CHKMEMQ;
1209a7e14dcfSSatish Balay   /*        [x1    (n)
1210a7e14dcfSSatish Balay              x2    (nb) may be 0
1211a7e14dcfSSatish Balay              x3    (me) may be 0
1212a7e14dcfSSatish Balay              x4    (nb) may be 0 */
1213a7e14dcfSSatish Balay   if (X1) {
1214a7e14dcfSSatish Balay     ierr = VecScatterBegin(ipmP->step1,STEP,X1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1215a7e14dcfSSatish Balay     ierr = VecScatterEnd(ipmP->step1,STEP,X1,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1216a7e14dcfSSatish Balay   }
1217a7e14dcfSSatish Balay   if (X2 && ipmP->nb > 0) {
1218a7e14dcfSSatish Balay     ierr = VecScatterBegin(ipmP->step2,STEP,X2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1219a7e14dcfSSatish Balay     ierr = VecScatterEnd(ipmP->step2,STEP,X2,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1220a7e14dcfSSatish Balay   }
1221a7e14dcfSSatish Balay   if (X3 && ipmP->me > 0) {
1222a7e14dcfSSatish Balay     ierr = VecScatterBegin(ipmP->step3,STEP,X3,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1223a7e14dcfSSatish Balay     ierr = VecScatterEnd(ipmP->step3,STEP,X3,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1224a7e14dcfSSatish Balay   }
1225a7e14dcfSSatish Balay   if (X4 && ipmP->nb > 0) {
1226a7e14dcfSSatish Balay     ierr = VecScatterBegin(ipmP->step4,STEP,X4,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1227a7e14dcfSSatish Balay     ierr = VecScatterEnd(ipmP->step4,STEP,X4,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1228a7e14dcfSSatish Balay   }
1229a7e14dcfSSatish Balay   CHKMEMQ;
1230a7e14dcfSSatish Balay #if defined(DEBUG_SCATTER)
1231a7e14dcfSSatish Balay   PetscPrintf(PETSC_COMM_WORLD,"Step\n");
1232a7e14dcfSSatish Balay   VecView(STEP,0);
1233a7e14dcfSSatish Balay   PetscPrintf(PETSC_COMM_WORLD,"X1-X4\n");
1234a7e14dcfSSatish Balay   if (X1) {VecView(X1,0);}
1235a7e14dcfSSatish Balay   if (X2) {VecView(X2,0);}
1236a7e14dcfSSatish Balay   if (X3) {VecView(X3,0);}
1237a7e14dcfSSatish Balay   if (X4) {VecView(X4,0);}
1238a7e14dcfSSatish Balay #endif
1239a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1240a7e14dcfSSatish Balay }
1241a7e14dcfSSatish Balay 
1242a7e14dcfSSatish Balay EXTERN_C_BEGIN
1243a7e14dcfSSatish Balay #undef __FUNCT__
1244a7e14dcfSSatish Balay #define __FUNCT__ "TaoCreate_IPM"
1245a7e14dcfSSatish Balay PetscErrorCode TaoCreate_IPM(TaoSolver tao)
1246a7e14dcfSSatish Balay {
1247a7e14dcfSSatish Balay   TAO_IPM        *ipmP;
1248a7e14dcfSSatish Balay   PetscErrorCode ierr;
1249a7e14dcfSSatish Balay 
1250a7e14dcfSSatish Balay   PetscFunctionBegin;
1251a7e14dcfSSatish Balay   tao->ops->setup = TaoSetup_IPM;
1252a7e14dcfSSatish Balay   tao->ops->solve = TaoSolve_IPM;
1253a7e14dcfSSatish Balay   tao->ops->view = TaoView_IPM;
1254a7e14dcfSSatish Balay   tao->ops->setfromoptions = TaoSetFromOptions_IPM;
1255a7e14dcfSSatish Balay   tao->ops->destroy = TaoDestroy_IPM;
1256e9f9aeaeSSatish Balay   /* tao->ops->computedual = TaoComputeDual_IPM; */
1257a7e14dcfSSatish Balay 
12583c9e27cfSGeoffrey Irving   ierr = PetscNewLog(tao,&ipmP);CHKERRQ(ierr);
1259a7e14dcfSSatish Balay   tao->data = (void*)ipmP;
1260a7e14dcfSSatish Balay   tao->max_it = 200;
1261a7e14dcfSSatish Balay   tao->max_funcs = 500;
1262a7e14dcfSSatish Balay   tao->fatol = 1e-4;
1263a7e14dcfSSatish Balay   tao->frtol = 1e-4;
1264a7e14dcfSSatish Balay   ipmP->dec = 10000; /* line search critera */
1265a7e14dcfSSatish Balay   ipmP->taumin = 0.995;
1266a7e14dcfSSatish Balay   ipmP->monitorkkt = PETSC_FALSE;
1267a7e14dcfSSatish Balay   ipmP->pushs = 100;
1268a7e14dcfSSatish Balay   ipmP->pushnu = 100;
1269a7e14dcfSSatish Balay   ierr = KSPCreate(((PetscObject)tao)->comm, &tao->ksp);CHKERRQ(ierr);
1270a7e14dcfSSatish Balay   PetscFunctionReturn(0);
1271a7e14dcfSSatish Balay }
1272a7e14dcfSSatish Balay EXTERN_C_END
1273a7e14dcfSSatish Balay 
1274