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