xref: /petsc/src/tao/constrained/impls/ipm/pdipm.c (revision aad13602431b2aae7f9a4915734a4ba0e9437f10)
1*aad13602SShrirang Abhyankar #include <petsctaolinesearch.h>
2*aad13602SShrirang Abhyankar #include <../src/tao/constrained/impls/ipm/pdipm.h>
3*aad13602SShrirang Abhyankar #include <petscsnes.h>
4*aad13602SShrirang Abhyankar 
5*aad13602SShrirang Abhyankar /*
6*aad13602SShrirang Abhyankar    TaoPDIPMEvaluateFunctionsAndJacobians - Evaluate the objective function f, gradient fx, constraints, and all the Jacobians at current vector
7*aad13602SShrirang Abhyankar 
8*aad13602SShrirang Abhyankar    Collective on tao
9*aad13602SShrirang Abhyankar 
10*aad13602SShrirang Abhyankar    Input Parameter:
11*aad13602SShrirang Abhyankar +  tao - solver context
12*aad13602SShrirang Abhyankar -  x - vector at which all objects to be evaluated
13*aad13602SShrirang Abhyankar 
14*aad13602SShrirang Abhyankar    Level: beginner
15*aad13602SShrirang Abhyankar 
16*aad13602SShrirang Abhyankar .seealso: TaoPDIPMUpdateConstraints(), TaoPDIPMSetUpBounds()
17*aad13602SShrirang Abhyankar */
18*aad13602SShrirang Abhyankar PetscErrorCode TaoPDIPMEvaluateFunctionsAndJacobians(Tao tao,Vec x)
19*aad13602SShrirang Abhyankar {
20*aad13602SShrirang Abhyankar   PetscErrorCode ierr;
21*aad13602SShrirang Abhyankar   TAO_PDIPM      *pdipm=(TAO_PDIPM*)tao->data;
22*aad13602SShrirang Abhyankar 
23*aad13602SShrirang Abhyankar   PetscFunctionBegin;
24*aad13602SShrirang Abhyankar   /* Compute user objective function and gradient */
25*aad13602SShrirang Abhyankar   ierr = TaoComputeObjectiveAndGradient(tao,x,&pdipm->obj,tao->gradient);CHKERRQ(ierr);
26*aad13602SShrirang Abhyankar 
27*aad13602SShrirang Abhyankar   /* Equality constraints and Jacobian */
28*aad13602SShrirang Abhyankar   if (pdipm->Ng) {
29*aad13602SShrirang Abhyankar     ierr = TaoComputeEqualityConstraints(tao,x,tao->constraints_equality);CHKERRQ(ierr);
30*aad13602SShrirang Abhyankar     ierr = TaoComputeJacobianEquality(tao,x,tao->jacobian_equality,tao->jacobian_equality_pre);CHKERRQ(ierr);
31*aad13602SShrirang Abhyankar   }
32*aad13602SShrirang Abhyankar 
33*aad13602SShrirang Abhyankar   /* Inequality constraints and Jacobian */
34*aad13602SShrirang Abhyankar   if (pdipm->Nh) {
35*aad13602SShrirang Abhyankar     ierr = TaoComputeInequalityConstraints(tao,x,tao->constraints_inequality);CHKERRQ(ierr);
36*aad13602SShrirang Abhyankar     ierr = TaoComputeJacobianInequality(tao,x,tao->jacobian_inequality,tao->jacobian_inequality_pre);CHKERRQ(ierr);
37*aad13602SShrirang Abhyankar   }
38*aad13602SShrirang Abhyankar   PetscFunctionReturn(0);
39*aad13602SShrirang Abhyankar }
40*aad13602SShrirang Abhyankar 
41*aad13602SShrirang Abhyankar /*
42*aad13602SShrirang Abhyankar   TaoPDIPMUpdateConstraints - Update the vectors ce and ci at x
43*aad13602SShrirang Abhyankar 
44*aad13602SShrirang Abhyankar   Collective
45*aad13602SShrirang Abhyankar 
46*aad13602SShrirang Abhyankar   Input Parameter:
47*aad13602SShrirang Abhyankar + tao - Tao context
48*aad13602SShrirang Abhyankar - x - vector at which constraints to be evaluted
49*aad13602SShrirang Abhyankar 
50*aad13602SShrirang Abhyankar    Level: beginner
51*aad13602SShrirang Abhyankar 
52*aad13602SShrirang Abhyankar .seealso: TaoPDIPMEvaluateFunctionsAndJacobians()
53*aad13602SShrirang Abhyankar */
54*aad13602SShrirang Abhyankar PetscErrorCode TaoPDIPMUpdateConstraints(Tao tao,Vec x)
55*aad13602SShrirang Abhyankar {
56*aad13602SShrirang Abhyankar   PetscErrorCode    ierr;
57*aad13602SShrirang Abhyankar   TAO_PDIPM         *pdipm=(TAO_PDIPM*)tao->data;
58*aad13602SShrirang Abhyankar   PetscInt          i,offset,offset1,k,xstart;
59*aad13602SShrirang Abhyankar   PetscScalar       *carr;
60*aad13602SShrirang Abhyankar   const PetscInt    *ubptr,*lbptr,*bxptr,*fxptr;
61*aad13602SShrirang Abhyankar   const PetscScalar *xarr,*xuarr,*xlarr,*garr,*harr;
62*aad13602SShrirang Abhyankar 
63*aad13602SShrirang Abhyankar   PetscFunctionBegin;
64*aad13602SShrirang Abhyankar   ierr = VecGetOwnershipRange(x,&xstart,NULL);CHKERRQ(ierr);
65*aad13602SShrirang Abhyankar 
66*aad13602SShrirang Abhyankar   ierr = VecGetArrayRead(x,&xarr);CHKERRQ(ierr);
67*aad13602SShrirang Abhyankar   ierr = VecGetArrayRead(tao->XU,&xuarr);CHKERRQ(ierr);
68*aad13602SShrirang Abhyankar   ierr = VecGetArrayRead(tao->XL,&xlarr);CHKERRQ(ierr);
69*aad13602SShrirang Abhyankar 
70*aad13602SShrirang Abhyankar   /* (1) Update ce vector */
71*aad13602SShrirang Abhyankar   ierr = VecGetArray(pdipm->ce,&carr);CHKERRQ(ierr);
72*aad13602SShrirang Abhyankar 
73*aad13602SShrirang Abhyankar   if(pdipm->Ng) {
74*aad13602SShrirang Abhyankar     /* (1.a) Inserting updated g(x) */
75*aad13602SShrirang Abhyankar     ierr = VecGetArrayRead(tao->constraints_equality,&garr);CHKERRQ(ierr);
76*aad13602SShrirang Abhyankar     ierr = PetscMemcpy(carr,garr,pdipm->ng*sizeof(PetscScalar));CHKERRQ(ierr);
77*aad13602SShrirang Abhyankar     ierr = VecRestoreArrayRead(tao->constraints_equality,&garr);CHKERRQ(ierr);
78*aad13602SShrirang Abhyankar   }
79*aad13602SShrirang Abhyankar 
80*aad13602SShrirang Abhyankar   /* (1.b) Update xfixed */
81*aad13602SShrirang Abhyankar   if (pdipm->Nxfixed) {
82*aad13602SShrirang Abhyankar     offset = pdipm->ng;
83*aad13602SShrirang Abhyankar     ierr = ISGetIndices(pdipm->isxfixed,&fxptr);CHKERRQ(ierr); /* global indices in x */
84*aad13602SShrirang Abhyankar     for (k=0;k < pdipm->nxfixed;k++){
85*aad13602SShrirang Abhyankar       i = fxptr[k]-xstart;
86*aad13602SShrirang Abhyankar       carr[offset + k] = xarr[i] - xuarr[i];
87*aad13602SShrirang Abhyankar     }
88*aad13602SShrirang Abhyankar   }
89*aad13602SShrirang Abhyankar   ierr = VecRestoreArray(pdipm->ce,&carr);CHKERRQ(ierr);
90*aad13602SShrirang Abhyankar 
91*aad13602SShrirang Abhyankar   /* (2) Update ci vector */
92*aad13602SShrirang Abhyankar   ierr = VecGetArray(pdipm->ci,&carr);CHKERRQ(ierr);
93*aad13602SShrirang Abhyankar 
94*aad13602SShrirang Abhyankar   if(pdipm->Nh) {
95*aad13602SShrirang Abhyankar     /* (2.a) Inserting updated h(x) */
96*aad13602SShrirang Abhyankar     ierr = VecGetArrayRead(tao->constraints_inequality,&harr);CHKERRQ(ierr);
97*aad13602SShrirang Abhyankar     ierr = PetscMemcpy(carr,harr,pdipm->nh*sizeof(PetscScalar));CHKERRQ(ierr);
98*aad13602SShrirang Abhyankar     ierr = VecRestoreArrayRead(tao->constraints_inequality,&harr);CHKERRQ(ierr);
99*aad13602SShrirang Abhyankar   }
100*aad13602SShrirang Abhyankar 
101*aad13602SShrirang Abhyankar   /* (2.b) Update xub */
102*aad13602SShrirang Abhyankar   offset = pdipm->nh;
103*aad13602SShrirang Abhyankar   if (pdipm->Nxub) {
104*aad13602SShrirang Abhyankar     ierr = ISGetIndices(pdipm->isxub,&ubptr);CHKERRQ(ierr);
105*aad13602SShrirang Abhyankar     for (k=0; k<pdipm->nxub; k++){
106*aad13602SShrirang Abhyankar       i = ubptr[k]-xstart;
107*aad13602SShrirang Abhyankar       carr[offset + k] = xuarr[i] - xarr[i];
108*aad13602SShrirang Abhyankar     }
109*aad13602SShrirang Abhyankar   }
110*aad13602SShrirang Abhyankar 
111*aad13602SShrirang Abhyankar   if (pdipm->Nxlb) {
112*aad13602SShrirang Abhyankar     /* (2.c) Update xlb */
113*aad13602SShrirang Abhyankar     offset += pdipm->nxub;
114*aad13602SShrirang Abhyankar     ierr = ISGetIndices(pdipm->isxlb,&lbptr);CHKERRQ(ierr); /* global indices in x */
115*aad13602SShrirang Abhyankar     for (k=0; k<pdipm->nxlb; k++){
116*aad13602SShrirang Abhyankar       i = lbptr[k]-xstart;
117*aad13602SShrirang Abhyankar       carr[offset + k] = xarr[i] - xlarr[i];
118*aad13602SShrirang Abhyankar     }
119*aad13602SShrirang Abhyankar   }
120*aad13602SShrirang Abhyankar 
121*aad13602SShrirang Abhyankar   if (pdipm->Nxbox) {
122*aad13602SShrirang Abhyankar     /* (2.d) Update xbox */
123*aad13602SShrirang Abhyankar     offset += pdipm->nxlb;
124*aad13602SShrirang Abhyankar     offset1 = offset + pdipm->nxbox;
125*aad13602SShrirang Abhyankar     ierr = ISGetIndices(pdipm->isxbox,&bxptr);CHKERRQ(ierr); /* global indices in x */
126*aad13602SShrirang Abhyankar     for (k=0; k<pdipm->nxbox; k++){
127*aad13602SShrirang Abhyankar       i = bxptr[k]-xstart; /* local indices in x */
128*aad13602SShrirang Abhyankar       carr[offset+k]  = xuarr[i] - xarr[i];
129*aad13602SShrirang Abhyankar       carr[offset1+k] = xarr[i]  - xlarr[i];
130*aad13602SShrirang Abhyankar     }
131*aad13602SShrirang Abhyankar   }
132*aad13602SShrirang Abhyankar   ierr = VecRestoreArray(pdipm->ci,&carr);CHKERRQ(ierr);
133*aad13602SShrirang Abhyankar 
134*aad13602SShrirang Abhyankar   /* Restoring Vectors */
135*aad13602SShrirang Abhyankar   ierr = VecRestoreArrayRead(x,&xarr);CHKERRQ(ierr);
136*aad13602SShrirang Abhyankar   ierr = VecRestoreArrayRead(tao->XU,&xuarr);CHKERRQ(ierr);
137*aad13602SShrirang Abhyankar   ierr = VecRestoreArrayRead(tao->XL,&xlarr);CHKERRQ(ierr);
138*aad13602SShrirang Abhyankar   PetscFunctionReturn(0);
139*aad13602SShrirang Abhyankar }
140*aad13602SShrirang Abhyankar 
141*aad13602SShrirang Abhyankar /*
142*aad13602SShrirang Abhyankar    TaoPDIPMSetUpBounds - Create upper and lower bound vectors of x
143*aad13602SShrirang Abhyankar 
144*aad13602SShrirang Abhyankar    Collective
145*aad13602SShrirang Abhyankar 
146*aad13602SShrirang Abhyankar    Input Parameter:
147*aad13602SShrirang Abhyankar .  tao - holds pdipm and XL & XU
148*aad13602SShrirang Abhyankar 
149*aad13602SShrirang Abhyankar    Level: beginner
150*aad13602SShrirang Abhyankar 
151*aad13602SShrirang Abhyankar .seealso: TaoPDIPMUpdateConstraints
152*aad13602SShrirang Abhyankar */
153*aad13602SShrirang Abhyankar PetscErrorCode TaoPDIPMSetUpBounds(Tao tao)
154*aad13602SShrirang Abhyankar {
155*aad13602SShrirang Abhyankar   PetscErrorCode    ierr;
156*aad13602SShrirang Abhyankar   TAO_PDIPM         *pdipm=(TAO_PDIPM*)tao->data;
157*aad13602SShrirang Abhyankar   const PetscScalar *xl,*xu;
158*aad13602SShrirang Abhyankar   PetscInt          n,*ixlb,*ixub,*ixfixed,*ixfree,*ixbox,i,low,high,idx;
159*aad13602SShrirang Abhyankar   MPI_Comm          comm;
160*aad13602SShrirang Abhyankar   PetscInt          sendbuf[5],recvbuf[5];
161*aad13602SShrirang Abhyankar 
162*aad13602SShrirang Abhyankar   PetscFunctionBegin;
163*aad13602SShrirang Abhyankar   /* Creates upper and lower bounds vectors on x, if not created already */
164*aad13602SShrirang Abhyankar   ierr = TaoComputeVariableBounds(tao);CHKERRQ(ierr);
165*aad13602SShrirang Abhyankar 
166*aad13602SShrirang Abhyankar   ierr = VecGetLocalSize(tao->XL,&n);CHKERRQ(ierr);
167*aad13602SShrirang Abhyankar   ierr = PetscMalloc5(n,&ixlb,n,&ixub,n,&ixfree,n,&ixfixed,n,&ixbox);CHKERRQ(ierr);
168*aad13602SShrirang Abhyankar 
169*aad13602SShrirang Abhyankar   ierr = VecGetOwnershipRange(tao->XL,&low,&high);CHKERRQ(ierr);
170*aad13602SShrirang Abhyankar   ierr = VecGetArrayRead(tao->XL,&xl);CHKERRQ(ierr);
171*aad13602SShrirang Abhyankar   ierr = VecGetArrayRead(tao->XU,&xu);CHKERRQ(ierr);
172*aad13602SShrirang Abhyankar   for (i=0; i<n; i++) {
173*aad13602SShrirang Abhyankar     idx = low + i;
174*aad13602SShrirang Abhyankar     if((PetscRealPart(xl[i]) > PETSC_NINFINITY) && (PetscRealPart(xu[i]) < PETSC_INFINITY)) {
175*aad13602SShrirang Abhyankar       if (PetscRealPart(xl[i]) == PetscRealPart(xu[i])) {
176*aad13602SShrirang Abhyankar         ixfixed[pdipm->nxfixed++]  = idx;
177*aad13602SShrirang Abhyankar       } else ixbox[pdipm->nxbox++] = idx;
178*aad13602SShrirang Abhyankar     } else {
179*aad13602SShrirang Abhyankar       if ((PetscRealPart(xl[i]) > PETSC_NINFINITY) && (PetscRealPart(xu[i]) >= PETSC_INFINITY)) {
180*aad13602SShrirang Abhyankar         ixlb[pdipm->nxlb++] = idx;
181*aad13602SShrirang Abhyankar       } else if ((PetscRealPart(xl[i]) <= PETSC_NINFINITY) && (PetscRealPart(xu[i]) < PETSC_INFINITY)) {
182*aad13602SShrirang Abhyankar         ixub[pdipm->nxlb++] = idx;
183*aad13602SShrirang Abhyankar       } else  ixfree[pdipm->nxfree++] = idx;
184*aad13602SShrirang Abhyankar     }
185*aad13602SShrirang Abhyankar   }
186*aad13602SShrirang Abhyankar   ierr = VecRestoreArrayRead(tao->XL,&xl);CHKERRQ(ierr);
187*aad13602SShrirang Abhyankar   ierr = VecRestoreArrayRead(tao->XU,&xu);CHKERRQ(ierr);
188*aad13602SShrirang Abhyankar 
189*aad13602SShrirang Abhyankar   ierr = PetscObjectGetComm((PetscObject)tao,&comm);CHKERRQ(ierr);
190*aad13602SShrirang Abhyankar   sendbuf[0] = pdipm->nxlb;
191*aad13602SShrirang Abhyankar   sendbuf[1] = pdipm->nxub;
192*aad13602SShrirang Abhyankar   sendbuf[2] = pdipm->nxfixed;
193*aad13602SShrirang Abhyankar   sendbuf[3] = pdipm->nxbox;
194*aad13602SShrirang Abhyankar   sendbuf[4] = pdipm->nxfree;
195*aad13602SShrirang Abhyankar 
196*aad13602SShrirang Abhyankar   ierr = MPI_Allreduce(sendbuf,recvbuf,5,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr);
197*aad13602SShrirang Abhyankar   pdipm->Nxlb    = recvbuf[0];
198*aad13602SShrirang Abhyankar   pdipm->Nxub    = recvbuf[1];
199*aad13602SShrirang Abhyankar   pdipm->Nxfixed = recvbuf[2];
200*aad13602SShrirang Abhyankar   pdipm->Nxbox   = recvbuf[3];
201*aad13602SShrirang Abhyankar   pdipm->Nxfree  = recvbuf[4];
202*aad13602SShrirang Abhyankar 
203*aad13602SShrirang Abhyankar   if (pdipm->Nxlb) {
204*aad13602SShrirang Abhyankar     ierr = ISCreateGeneral(comm,pdipm->nxlb,ixlb,PETSC_COPY_VALUES,&pdipm->isxlb);CHKERRQ(ierr);
205*aad13602SShrirang Abhyankar   }
206*aad13602SShrirang Abhyankar   if (pdipm->Nxub) {
207*aad13602SShrirang Abhyankar     ierr = ISCreateGeneral(comm,pdipm->nxub,ixub,PETSC_COPY_VALUES,&pdipm->isxub);CHKERRQ(ierr);
208*aad13602SShrirang Abhyankar   }
209*aad13602SShrirang Abhyankar   if (pdipm->Nxfixed) {
210*aad13602SShrirang Abhyankar     ierr = ISCreateGeneral(comm,pdipm->nxfixed,ixfixed,PETSC_COPY_VALUES,&pdipm->isxfixed);CHKERRQ(ierr);
211*aad13602SShrirang Abhyankar   }
212*aad13602SShrirang Abhyankar   if (pdipm->Nxbox) {
213*aad13602SShrirang Abhyankar     ierr = ISCreateGeneral(comm,pdipm->nxbox,ixbox,PETSC_COPY_VALUES,&pdipm->isxbox);CHKERRQ(ierr);
214*aad13602SShrirang Abhyankar   }
215*aad13602SShrirang Abhyankar   if (pdipm->Nxfree) {
216*aad13602SShrirang Abhyankar     ierr = ISCreateGeneral(comm,pdipm->nxfree,ixfree,PETSC_COPY_VALUES,&pdipm->isxfree);CHKERRQ(ierr);
217*aad13602SShrirang Abhyankar   }
218*aad13602SShrirang Abhyankar   ierr = PetscFree5(ixlb,ixub,ixfixed,ixbox,ixfree);CHKERRQ(ierr);
219*aad13602SShrirang Abhyankar   PetscFunctionReturn(0);
220*aad13602SShrirang Abhyankar }
221*aad13602SShrirang Abhyankar 
222*aad13602SShrirang Abhyankar /*
223*aad13602SShrirang Abhyankar    TaoPDIPMInitializeSolution - Initialize PDIPM solution X = [x; lambdae; lambdai; z].
224*aad13602SShrirang Abhyankar    X consists of four subvectors in the order [x; lambdae; lambdai; z]. These
225*aad13602SShrirang Abhyankar      four subvectors need to be initialized and its values copied over to X. Instead
226*aad13602SShrirang Abhyankar      of copying, we use VecPlace/ResetArray functions to share the memory locations for
227*aad13602SShrirang Abhyankar      X and the subvectors
228*aad13602SShrirang Abhyankar 
229*aad13602SShrirang Abhyankar    Collective
230*aad13602SShrirang Abhyankar 
231*aad13602SShrirang Abhyankar    Input Parameter:
232*aad13602SShrirang Abhyankar .  tao - Tao context
233*aad13602SShrirang Abhyankar 
234*aad13602SShrirang Abhyankar    Level: beginner
235*aad13602SShrirang Abhyankar */
236*aad13602SShrirang Abhyankar PetscErrorCode TaoPDIPMInitializeSolution(Tao tao)
237*aad13602SShrirang Abhyankar {
238*aad13602SShrirang Abhyankar   PetscErrorCode ierr;
239*aad13602SShrirang Abhyankar   TAO_PDIPM      *pdipm = (TAO_PDIPM*)tao->data;
240*aad13602SShrirang Abhyankar   PetscScalar    *Xarr,*z,*lambdai;
241*aad13602SShrirang Abhyankar   PetscInt       i;
242*aad13602SShrirang Abhyankar   const PetscScalar *xarr,*h;
243*aad13602SShrirang Abhyankar 
244*aad13602SShrirang Abhyankar   PetscFunctionBegin;
245*aad13602SShrirang Abhyankar   ierr = VecGetArray(pdipm->X,&Xarr);CHKERRQ(ierr);
246*aad13602SShrirang Abhyankar 
247*aad13602SShrirang Abhyankar   /* Set Initialize X.x = tao->solution */
248*aad13602SShrirang Abhyankar   ierr = VecGetArrayRead(tao->solution,&xarr);CHKERRQ(ierr);
249*aad13602SShrirang Abhyankar   ierr = PetscMemcpy(Xarr,xarr,pdipm->nx*sizeof(PetscScalar));CHKERRQ(ierr);
250*aad13602SShrirang Abhyankar   ierr = VecRestoreArrayRead(tao->solution,&xarr);CHKERRQ(ierr);
251*aad13602SShrirang Abhyankar 
252*aad13602SShrirang Abhyankar   /* Initialize X.lambdae = 0.0 */
253*aad13602SShrirang Abhyankar   ierr = VecSet(pdipm->lambdae,0.0);CHKERRQ(ierr);
254*aad13602SShrirang Abhyankar 
255*aad13602SShrirang Abhyankar   /* Initialize X.lambdai = push_init_lambdai, X.z = push_init_slack */
256*aad13602SShrirang Abhyankar   ierr = VecSet(pdipm->lambdai,pdipm->push_init_lambdai);CHKERRQ(ierr);
257*aad13602SShrirang Abhyankar   ierr = VecSet(pdipm->z,pdipm->push_init_slack);CHKERRQ(ierr);
258*aad13602SShrirang Abhyankar 
259*aad13602SShrirang Abhyankar   /* Additional modification for X.lambdai and X.z */
260*aad13602SShrirang Abhyankar   ierr = VecGetArray(pdipm->lambdai,&lambdai);CHKERRQ(ierr);
261*aad13602SShrirang Abhyankar   ierr = VecGetArray(pdipm->z,&z);CHKERRQ(ierr);
262*aad13602SShrirang Abhyankar   if(pdipm->Nh) {
263*aad13602SShrirang Abhyankar     ierr = VecGetArrayRead(tao->constraints_inequality,&h);CHKERRQ(ierr);
264*aad13602SShrirang Abhyankar     for (i=0; i < pdipm->nh; i++) {
265*aad13602SShrirang Abhyankar       if (h[i] < -pdipm->push_init_slack) z[i] = -h[i];
266*aad13602SShrirang Abhyankar       if (pdipm->mu/z[i] > pdipm->push_init_lambdai) lambdai[i] = pdipm->mu/z[i];
267*aad13602SShrirang Abhyankar     }
268*aad13602SShrirang Abhyankar     ierr = VecRestoreArrayRead(tao->constraints_inequality,&h);CHKERRQ(ierr);
269*aad13602SShrirang Abhyankar   }
270*aad13602SShrirang Abhyankar   ierr = VecRestoreArray(pdipm->lambdai,&lambdai);CHKERRQ(ierr);
271*aad13602SShrirang Abhyankar   ierr = VecRestoreArray(pdipm->z,&z);CHKERRQ(ierr);
272*aad13602SShrirang Abhyankar 
273*aad13602SShrirang Abhyankar   ierr = VecRestoreArray(pdipm->X,&Xarr);CHKERRQ(ierr);
274*aad13602SShrirang Abhyankar   PetscFunctionReturn(0);
275*aad13602SShrirang Abhyankar }
276*aad13602SShrirang Abhyankar 
277*aad13602SShrirang Abhyankar /*
278*aad13602SShrirang Abhyankar    TaoSNESJacobian_PDIPM - Evaluate the Hessian matrix at X
279*aad13602SShrirang Abhyankar 
280*aad13602SShrirang Abhyankar    Input Parameter:
281*aad13602SShrirang Abhyankar    snes - SNES context
282*aad13602SShrirang Abhyankar    X - KKT Vector
283*aad13602SShrirang Abhyankar    *ctx - pdipm context
284*aad13602SShrirang Abhyankar 
285*aad13602SShrirang Abhyankar    Output Parameter:
286*aad13602SShrirang Abhyankar    J - Hessian matrix
287*aad13602SShrirang Abhyankar    Jpre - Preconditioner
288*aad13602SShrirang Abhyankar */
289*aad13602SShrirang Abhyankar PetscErrorCode TaoSNESJacobian_PDIPM(SNES snes,Vec X, Mat J, Mat Jpre, void *ctx)
290*aad13602SShrirang Abhyankar {
291*aad13602SShrirang Abhyankar   PetscErrorCode    ierr;
292*aad13602SShrirang Abhyankar   Tao               tao=(Tao)ctx;
293*aad13602SShrirang Abhyankar   TAO_PDIPM         *pdipm = (TAO_PDIPM*)tao->data;
294*aad13602SShrirang Abhyankar   PetscInt          i,row,cols[2],Jrstart,rjstart,nc,j;
295*aad13602SShrirang Abhyankar   const PetscInt    *aj,*ranges,*Jranges,*rranges,*cranges;
296*aad13602SShrirang Abhyankar   const PetscScalar *Xarr,*aa;
297*aad13602SShrirang Abhyankar   PetscScalar       vals[2];
298*aad13602SShrirang Abhyankar   PetscInt          proc,nx_all,*nce_all=pdipm->nce_all;
299*aad13602SShrirang Abhyankar   MPI_Comm          comm;
300*aad13602SShrirang Abhyankar   PetscMPIInt       rank,size;
301*aad13602SShrirang Abhyankar   Mat               jac_equality_trans=pdipm->jac_equality_trans,jac_inequality_trans=pdipm->jac_inequality_trans;
302*aad13602SShrirang Abhyankar 
303*aad13602SShrirang Abhyankar   PetscFunctionBegin;
304*aad13602SShrirang Abhyankar   ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
305*aad13602SShrirang Abhyankar   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
306*aad13602SShrirang Abhyankar   ierr = MPI_Comm_rank(comm,&size);CHKERRQ(ierr);
307*aad13602SShrirang Abhyankar 
308*aad13602SShrirang Abhyankar   ierr = MatGetOwnershipRanges(Jpre,&Jranges);CHKERRQ(ierr);
309*aad13602SShrirang Abhyankar   ierr = MatGetOwnershipRange(Jpre,&Jrstart,NULL);CHKERRQ(ierr);
310*aad13602SShrirang Abhyankar   ierr = MatGetOwnershipRangesColumn(tao->hessian,&rranges);CHKERRQ(ierr);
311*aad13602SShrirang Abhyankar   ierr = MatGetOwnershipRangesColumn(tao->hessian,&cranges);CHKERRQ(ierr);
312*aad13602SShrirang Abhyankar 
313*aad13602SShrirang Abhyankar   ierr = VecGetArrayRead(X,&Xarr);CHKERRQ(ierr);
314*aad13602SShrirang Abhyankar 
315*aad13602SShrirang Abhyankar   /* (2) insert Z and Ci to Jpre -- overwrite existing values */
316*aad13602SShrirang Abhyankar   for (i=0; i < pdipm->nci; i++) {
317*aad13602SShrirang Abhyankar     row     = Jrstart + pdipm->off_z + i;
318*aad13602SShrirang Abhyankar     cols[0] = Jrstart + pdipm->off_lambdai + i;
319*aad13602SShrirang Abhyankar     cols[1] = row;
320*aad13602SShrirang Abhyankar     vals[0] = Xarr[pdipm->off_z + i];
321*aad13602SShrirang Abhyankar     vals[1] = Xarr[pdipm->off_lambdai + i];
322*aad13602SShrirang Abhyankar     ierr = MatSetValues(Jpre,1,&row,2,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
323*aad13602SShrirang Abhyankar   }
324*aad13602SShrirang Abhyankar 
325*aad13602SShrirang Abhyankar   /* (3) insert 2nd row block of Jpre: [ grad g, 0, 0, 0] */
326*aad13602SShrirang Abhyankar   if(pdipm->Ng) {
327*aad13602SShrirang Abhyankar     ierr = MatGetOwnershipRange(tao->jacobian_equality,&rjstart,NULL);CHKERRQ(ierr);
328*aad13602SShrirang Abhyankar     for (i=0; i<pdipm->ng; i++){
329*aad13602SShrirang Abhyankar       row = Jrstart + pdipm->off_lambdae + i;
330*aad13602SShrirang Abhyankar 
331*aad13602SShrirang Abhyankar       ierr = MatGetRow(tao->jacobian_equality,i+rjstart,&nc,&aj,&aa);CHKERRQ(ierr);
332*aad13602SShrirang Abhyankar       proc = 0;
333*aad13602SShrirang Abhyankar       for (j=0; j < nc; j++) {
334*aad13602SShrirang Abhyankar         while (aj[j] >= cranges[proc+1]) proc++;
335*aad13602SShrirang Abhyankar         cols[0] = aj[j] - cranges[proc] + Jranges[proc];
336*aad13602SShrirang Abhyankar         ierr = MatSetValue(Jpre,row,cols[0],aa[j],INSERT_VALUES);CHKERRQ(ierr);
337*aad13602SShrirang Abhyankar       }
338*aad13602SShrirang Abhyankar       ierr = MatRestoreRow(tao->jacobian_equality,i+rjstart,&nc,&aj,&aa);CHKERRQ(ierr);
339*aad13602SShrirang Abhyankar     }
340*aad13602SShrirang Abhyankar   }
341*aad13602SShrirang Abhyankar 
342*aad13602SShrirang Abhyankar   if(pdipm->Nh) {
343*aad13602SShrirang Abhyankar     /* (4) insert 3nd row block of Jpre: [ grad h, 0, 0, 0] */
344*aad13602SShrirang Abhyankar     ierr = MatGetOwnershipRange(tao->jacobian_inequality,&rjstart,NULL);CHKERRQ(ierr);
345*aad13602SShrirang Abhyankar     for (i=0; i < pdipm->nh; i++){
346*aad13602SShrirang Abhyankar       row = Jrstart + pdipm->off_lambdai + i;
347*aad13602SShrirang Abhyankar 
348*aad13602SShrirang Abhyankar       ierr = MatGetRow(tao->jacobian_inequality,i+rjstart,&nc,&aj,&aa);CHKERRQ(ierr);
349*aad13602SShrirang Abhyankar       proc = 0;
350*aad13602SShrirang Abhyankar       for (j=0; j < nc; j++) {
351*aad13602SShrirang Abhyankar         while (aj[j] >= cranges[proc+1]) proc++;
352*aad13602SShrirang Abhyankar         cols[0] = aj[j] - cranges[proc] + Jranges[proc];
353*aad13602SShrirang Abhyankar         ierr = MatSetValue(Jpre,row,cols[0],aa[j],INSERT_VALUES);CHKERRQ(ierr);
354*aad13602SShrirang Abhyankar       }
355*aad13602SShrirang Abhyankar     ierr = MatRestoreRow(tao->jacobian_inequality,i+rjstart,&nc,&aj,&aa);CHKERRQ(ierr);
356*aad13602SShrirang Abhyankar     }
357*aad13602SShrirang Abhyankar   }
358*aad13602SShrirang Abhyankar 
359*aad13602SShrirang Abhyankar   /* (5) insert Wxx, grad g' and -grad h' to Jpre */
360*aad13602SShrirang Abhyankar   if(pdipm->Ng) {
361*aad13602SShrirang Abhyankar     ierr = MatTranspose(tao->jacobian_equality,MAT_REUSE_MATRIX,&jac_equality_trans);CHKERRQ(ierr);
362*aad13602SShrirang Abhyankar   }
363*aad13602SShrirang Abhyankar   if(pdipm->Nh) {
364*aad13602SShrirang Abhyankar     ierr = MatTranspose(tao->jacobian_inequality,MAT_REUSE_MATRIX,&jac_inequality_trans);CHKERRQ(ierr);
365*aad13602SShrirang Abhyankar   }
366*aad13602SShrirang Abhyankar 
367*aad13602SShrirang Abhyankar   ierr = VecPlaceArray(pdipm->x,Xarr);CHKERRQ(ierr);
368*aad13602SShrirang Abhyankar   ierr = TaoComputeHessian(tao,pdipm->x,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
369*aad13602SShrirang Abhyankar   ierr = VecResetArray(pdipm->x);CHKERRQ(ierr);
370*aad13602SShrirang Abhyankar 
371*aad13602SShrirang Abhyankar   ierr = MatGetOwnershipRange(tao->hessian,&rjstart,NULL);CHKERRQ(ierr);
372*aad13602SShrirang Abhyankar   for (i=0; i<pdipm->nx; i++){
373*aad13602SShrirang Abhyankar     row = Jrstart + i;
374*aad13602SShrirang Abhyankar 
375*aad13602SShrirang Abhyankar     /* insert Wxx */
376*aad13602SShrirang Abhyankar     ierr = MatGetRow(tao->hessian,i+rjstart,&nc,&aj,&aa);CHKERRQ(ierr);
377*aad13602SShrirang Abhyankar     proc = 0;
378*aad13602SShrirang Abhyankar     for (j=0; j < nc; j++) {
379*aad13602SShrirang Abhyankar       while (aj[j] >= cranges[proc+1]) proc++;
380*aad13602SShrirang Abhyankar       cols[0] = aj[j] - cranges[proc] + Jranges[proc];
381*aad13602SShrirang Abhyankar       ierr = MatSetValue(Jpre,row,cols[0],aa[j],INSERT_VALUES);CHKERRQ(ierr);
382*aad13602SShrirang Abhyankar     }
383*aad13602SShrirang Abhyankar     ierr = MatRestoreRow(tao->hessian,i+rjstart,&nc,&aj,&aa);CHKERRQ(ierr);
384*aad13602SShrirang Abhyankar 
385*aad13602SShrirang Abhyankar     if(pdipm->ng) {
386*aad13602SShrirang Abhyankar       /* insert grad g' */
387*aad13602SShrirang Abhyankar       ierr = MatGetRow(jac_equality_trans,i+rjstart,&nc,&aj,&aa);CHKERRQ(ierr);
388*aad13602SShrirang Abhyankar       ierr = MatGetOwnershipRanges(tao->jacobian_equality,&ranges);CHKERRQ(ierr);
389*aad13602SShrirang Abhyankar       proc = 0;
390*aad13602SShrirang Abhyankar       for (j=0; j < nc; j++) {
391*aad13602SShrirang Abhyankar         /* find row ownership of */
392*aad13602SShrirang Abhyankar         while (aj[j] >= ranges[proc+1]) proc++;
393*aad13602SShrirang Abhyankar         nx_all = rranges[proc+1] - rranges[proc];
394*aad13602SShrirang Abhyankar         cols[0] = aj[j] - ranges[proc] + Jranges[proc] + nx_all;
395*aad13602SShrirang Abhyankar         ierr = MatSetValue(Jpre,row,cols[0],aa[j],INSERT_VALUES);CHKERRQ(ierr);
396*aad13602SShrirang Abhyankar       }
397*aad13602SShrirang Abhyankar       ierr = MatRestoreRow(jac_equality_trans,i+rjstart,&nc,&aj,&aa);CHKERRQ(ierr);
398*aad13602SShrirang Abhyankar     }
399*aad13602SShrirang Abhyankar 
400*aad13602SShrirang Abhyankar     if(pdipm->nh) {
401*aad13602SShrirang Abhyankar       /* insert -grad h' */
402*aad13602SShrirang Abhyankar       ierr = MatGetRow(jac_inequality_trans,i+rjstart,&nc,&aj,&aa);CHKERRQ(ierr);
403*aad13602SShrirang Abhyankar       ierr = MatGetOwnershipRanges(tao->jacobian_inequality,&ranges);CHKERRQ(ierr);
404*aad13602SShrirang Abhyankar       proc = 0;
405*aad13602SShrirang Abhyankar       for (j=0; j < nc; j++) {
406*aad13602SShrirang Abhyankar         /* find row ownership of */
407*aad13602SShrirang Abhyankar         while (aj[j] >= ranges[proc+1]) proc++;
408*aad13602SShrirang Abhyankar         nx_all = rranges[proc+1] - rranges[proc];
409*aad13602SShrirang Abhyankar         cols[0] = aj[j] - ranges[proc] + Jranges[proc] + nx_all + nce_all[proc];
410*aad13602SShrirang Abhyankar         ierr = MatSetValue(Jpre,row,cols[0],-aa[j],INSERT_VALUES);CHKERRQ(ierr);
411*aad13602SShrirang Abhyankar       }
412*aad13602SShrirang Abhyankar       ierr = MatRestoreRow(jac_inequality_trans,i+rjstart,&nc,&aj,&aa);CHKERRQ(ierr);
413*aad13602SShrirang Abhyankar     }
414*aad13602SShrirang Abhyankar   }
415*aad13602SShrirang Abhyankar   ierr = VecRestoreArrayRead(X,&Xarr);CHKERRQ(ierr);
416*aad13602SShrirang Abhyankar 
417*aad13602SShrirang Abhyankar   /* (6) assemble Jpre and J */
418*aad13602SShrirang Abhyankar   ierr = MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
419*aad13602SShrirang Abhyankar   ierr = MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
420*aad13602SShrirang Abhyankar 
421*aad13602SShrirang Abhyankar   if (J != Jpre) {
422*aad13602SShrirang Abhyankar     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
423*aad13602SShrirang Abhyankar     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
424*aad13602SShrirang Abhyankar   }
425*aad13602SShrirang Abhyankar   PetscFunctionReturn(0);
426*aad13602SShrirang Abhyankar }
427*aad13602SShrirang Abhyankar 
428*aad13602SShrirang Abhyankar /*
429*aad13602SShrirang Abhyankar    TaoSnesFunction_PDIPM - Evaluate KKT function at X
430*aad13602SShrirang Abhyankar 
431*aad13602SShrirang Abhyankar    Input Parameter:
432*aad13602SShrirang Abhyankar    snes - SNES context
433*aad13602SShrirang Abhyankar    X - KKT Vector
434*aad13602SShrirang Abhyankar    *ctx - pdipm
435*aad13602SShrirang Abhyankar 
436*aad13602SShrirang Abhyankar    Output Parameter:
437*aad13602SShrirang Abhyankar    F - Updated Lagrangian vector
438*aad13602SShrirang Abhyankar */
439*aad13602SShrirang Abhyankar PetscErrorCode TaoSNESFunction_PDIPM(SNES snes,Vec X,Vec F,void *ctx)
440*aad13602SShrirang Abhyankar {
441*aad13602SShrirang Abhyankar   PetscErrorCode ierr;
442*aad13602SShrirang Abhyankar   Tao            tao=(Tao)ctx;
443*aad13602SShrirang Abhyankar   TAO_PDIPM      *pdipm = (TAO_PDIPM*)tao->data;
444*aad13602SShrirang Abhyankar   PetscScalar    *Farr;
445*aad13602SShrirang Abhyankar   Vec            x,L1;
446*aad13602SShrirang Abhyankar   PetscInt       i;
447*aad13602SShrirang Abhyankar   PetscReal      res[2],cnorm[2];
448*aad13602SShrirang Abhyankar   const PetscScalar *Xarr,*carr,*zarr,*larr;
449*aad13602SShrirang Abhyankar 
450*aad13602SShrirang Abhyankar   PetscFunctionBegin;
451*aad13602SShrirang Abhyankar   ierr = VecSet(F,0.0);CHKERRQ(ierr);
452*aad13602SShrirang Abhyankar 
453*aad13602SShrirang Abhyankar   ierr = VecGetArrayRead(X,&Xarr);CHKERRQ(ierr);
454*aad13602SShrirang Abhyankar   ierr = VecGetArray(F,&Farr);CHKERRQ(ierr);
455*aad13602SShrirang Abhyankar 
456*aad13602SShrirang Abhyankar   /* (0) Evaluate f, fx, Gx, Hx at X.x Note: pdipm->x is not changed below */
457*aad13602SShrirang Abhyankar   x = pdipm->x;
458*aad13602SShrirang Abhyankar   ierr = VecPlaceArray(x,Xarr);CHKERRQ(ierr);
459*aad13602SShrirang Abhyankar   ierr = TaoPDIPMEvaluateFunctionsAndJacobians(tao,x);CHKERRQ(ierr);
460*aad13602SShrirang Abhyankar 
461*aad13602SShrirang Abhyankar   /* Update ce, ci, and Jci at X.x */
462*aad13602SShrirang Abhyankar   ierr = TaoPDIPMUpdateConstraints(tao,x);CHKERRQ(ierr);
463*aad13602SShrirang Abhyankar   ierr = VecResetArray(x);CHKERRQ(ierr);
464*aad13602SShrirang Abhyankar 
465*aad13602SShrirang Abhyankar   /* (1) L1 = fx + (gradG'*DE + Jce_xfixed'*lambdae_xfixed) - (gradH'*DI + Jci_xb'*lambdai_xb) */
466*aad13602SShrirang Abhyankar   L1 = pdipm->x;
467*aad13602SShrirang Abhyankar   ierr = VecPlaceArray(L1,Farr);CHKERRQ(ierr);
468*aad13602SShrirang Abhyankar   if (pdipm->Nci) {
469*aad13602SShrirang Abhyankar     if(pdipm->Nh) {
470*aad13602SShrirang Abhyankar       /* L1 += gradH'*DI. Note: tao->DI is not changed below */
471*aad13602SShrirang Abhyankar       ierr = VecPlaceArray(tao->DI,Xarr+pdipm->off_lambdai);CHKERRQ(ierr);
472*aad13602SShrirang Abhyankar       ierr = MatMultTransposeAdd(tao->jacobian_inequality,tao->DI,L1,L1);CHKERRQ(ierr);
473*aad13602SShrirang Abhyankar       ierr = VecResetArray(tao->DI);CHKERRQ(ierr);
474*aad13602SShrirang Abhyankar     }
475*aad13602SShrirang Abhyankar 
476*aad13602SShrirang Abhyankar     /* L1 += Jci_xb'*lambdai_xb */
477*aad13602SShrirang Abhyankar     ierr = VecPlaceArray(pdipm->lambdai_xb,Xarr+pdipm->off_lambdai+pdipm->nh);CHKERRQ(ierr);
478*aad13602SShrirang Abhyankar     ierr = MatMultTransposeAdd(pdipm->Jci_xb,pdipm->lambdai_xb,L1,L1);CHKERRQ(ierr);
479*aad13602SShrirang Abhyankar     ierr = VecResetArray(pdipm->lambdai_xb);CHKERRQ(ierr);
480*aad13602SShrirang Abhyankar 
481*aad13602SShrirang Abhyankar     /* (1.4) L1 = - (gradH'*DI + Jci_xb'*lambdai_xb) */
482*aad13602SShrirang Abhyankar     ierr = VecScale(L1,-1.0);CHKERRQ(ierr);
483*aad13602SShrirang Abhyankar   }
484*aad13602SShrirang Abhyankar 
485*aad13602SShrirang Abhyankar   /* L1 += fx */
486*aad13602SShrirang Abhyankar   ierr = VecAXPY(L1,1.0,tao->gradient);CHKERRQ(ierr);
487*aad13602SShrirang Abhyankar 
488*aad13602SShrirang Abhyankar   if (pdipm->Nce) {
489*aad13602SShrirang Abhyankar     if(pdipm->Ng) {
490*aad13602SShrirang Abhyankar       /* L1 += gradG'*DE. Note: tao->DE is not changed below */
491*aad13602SShrirang Abhyankar       ierr = VecPlaceArray(tao->DE,Xarr+pdipm->off_lambdae);CHKERRQ(ierr);
492*aad13602SShrirang Abhyankar       ierr = MatMultTransposeAdd(tao->jacobian_equality,tao->DE,L1,L1);CHKERRQ(ierr);
493*aad13602SShrirang Abhyankar       ierr = VecResetArray(tao->DE);CHKERRQ(ierr);
494*aad13602SShrirang Abhyankar     }
495*aad13602SShrirang Abhyankar     if (pdipm->Nxfixed) {
496*aad13602SShrirang Abhyankar       /* L1 += Jce_xfixed'*lambdae_xfixed */
497*aad13602SShrirang Abhyankar       ierr = VecPlaceArray(pdipm->lambdae_xfixed,Xarr+pdipm->off_lambdae+pdipm->ng);CHKERRQ(ierr);
498*aad13602SShrirang Abhyankar       ierr = MatMultTransposeAdd(pdipm->Jce_xfixed,pdipm->lambdae_xfixed,L1,L1);CHKERRQ(ierr);
499*aad13602SShrirang Abhyankar       ierr = VecResetArray(pdipm->lambdae_xfixed);CHKERRQ(ierr);
500*aad13602SShrirang Abhyankar     }
501*aad13602SShrirang Abhyankar   }
502*aad13602SShrirang Abhyankar   ierr = VecNorm(L1,NORM_2,&res[0]);CHKERRQ(ierr);
503*aad13602SShrirang Abhyankar   ierr = VecResetArray(L1);CHKERRQ(ierr);
504*aad13602SShrirang Abhyankar 
505*aad13602SShrirang Abhyankar   /* (2) L2 = ce(x) */
506*aad13602SShrirang Abhyankar   if (pdipm->Nce) {
507*aad13602SShrirang Abhyankar     ierr = VecGetArrayRead(pdipm->ce,&carr);CHKERRQ(ierr);
508*aad13602SShrirang Abhyankar     for (i=0; i<pdipm->nce; i++) Farr[pdipm->off_lambdae + i] = carr[i];
509*aad13602SShrirang Abhyankar     ierr = VecRestoreArrayRead(pdipm->ce,&carr);CHKERRQ(ierr);
510*aad13602SShrirang Abhyankar   }
511*aad13602SShrirang Abhyankar   ierr = VecNorm(pdipm->ce,NORM_2,&cnorm[0]);CHKERRQ(ierr);
512*aad13602SShrirang Abhyankar 
513*aad13602SShrirang Abhyankar   if (pdipm->Nci) {
514*aad13602SShrirang Abhyankar     /* (3) L3 = ci(x) - z;
515*aad13602SShrirang Abhyankar        (4) L4 = Z * Lambdai * e - mu * e
516*aad13602SShrirang Abhyankar     */
517*aad13602SShrirang Abhyankar     ierr = VecGetArrayRead(pdipm->ci,&carr);CHKERRQ(ierr);
518*aad13602SShrirang Abhyankar     larr = Xarr+pdipm->off_lambdai;
519*aad13602SShrirang Abhyankar     zarr = Xarr+pdipm->off_z;
520*aad13602SShrirang Abhyankar     for (i=0; i<pdipm->nci; i++) {
521*aad13602SShrirang Abhyankar       Farr[pdipm->off_lambdai + i] = carr[i] - zarr[i];
522*aad13602SShrirang Abhyankar       Farr[pdipm->off_z       + i] = zarr[i]*larr[i] - pdipm->mu;
523*aad13602SShrirang Abhyankar     }
524*aad13602SShrirang Abhyankar     ierr = VecRestoreArrayRead(pdipm->ci,&carr);CHKERRQ(ierr);
525*aad13602SShrirang Abhyankar   }
526*aad13602SShrirang Abhyankar 
527*aad13602SShrirang Abhyankar   ierr = VecPlaceArray(pdipm->ci,Farr+pdipm->off_lambdai);CHKERRQ(ierr);
528*aad13602SShrirang Abhyankar   ierr = VecNorm(pdipm->ci,NORM_2,&cnorm[1]);CHKERRQ(ierr);
529*aad13602SShrirang Abhyankar   ierr = VecResetArray(pdipm->ci);CHKERRQ(ierr);
530*aad13602SShrirang Abhyankar 
531*aad13602SShrirang Abhyankar   /* note: pdipm->z is not changed below */
532*aad13602SShrirang Abhyankar   ierr = VecPlaceArray(pdipm->z,Farr+pdipm->off_z);CHKERRQ(ierr);
533*aad13602SShrirang Abhyankar   ierr = VecNorm(pdipm->z,NORM_2,&res[1]);CHKERRQ(ierr);
534*aad13602SShrirang Abhyankar   ierr = VecResetArray(pdipm->z);CHKERRQ(ierr);
535*aad13602SShrirang Abhyankar 
536*aad13602SShrirang Abhyankar   tao->residual = PetscSqrtReal(res[0]*res[0] + res[1]*res[1]);
537*aad13602SShrirang Abhyankar   tao->cnorm    = PetscSqrtReal(cnorm[0]*cnorm[0] + cnorm[1]*cnorm[1]);
538*aad13602SShrirang Abhyankar 
539*aad13602SShrirang Abhyankar   ierr = VecRestoreArrayRead(X,&Xarr);CHKERRQ(ierr);
540*aad13602SShrirang Abhyankar   ierr = VecRestoreArray(F,&Farr);CHKERRQ(ierr);
541*aad13602SShrirang Abhyankar   PetscFunctionReturn(0);
542*aad13602SShrirang Abhyankar }
543*aad13602SShrirang Abhyankar 
544*aad13602SShrirang Abhyankar /*
545*aad13602SShrirang Abhyankar    PDIPMLineSearch - Custom line search used with PDIPM.
546*aad13602SShrirang Abhyankar 
547*aad13602SShrirang Abhyankar    Collective on TAO
548*aad13602SShrirang Abhyankar 
549*aad13602SShrirang Abhyankar    Notes:
550*aad13602SShrirang Abhyankar    PDIPMLineSearch employs a simple backtracking line-search to keep
551*aad13602SShrirang Abhyankar    the slack variables (z) and inequality constraints lagrange multipliers
552*aad13602SShrirang Abhyankar    (lambdai) positive, i.e., z,lambdai >=0. It does this by calculating scalars
553*aad13602SShrirang Abhyankar    alpha_p and alpha_d to keep z,lambdai non-negative. The decision (x), and the
554*aad13602SShrirang Abhyankar    slack variables are updated as X = X + alpha_d*dx. The constraint multipliers
555*aad13602SShrirang Abhyankar    are updated as Lambdai = Lambdai + alpha_p*dLambdai. The barrier parameter mu
556*aad13602SShrirang Abhyankar    is also updated as mu = mu + z'lambdai/Nci
557*aad13602SShrirang Abhyankar */
558*aad13602SShrirang Abhyankar PetscErrorCode PDIPMLineSearch(SNESLineSearch linesearch,void *ctx)
559*aad13602SShrirang Abhyankar {
560*aad13602SShrirang Abhyankar   PetscErrorCode ierr;
561*aad13602SShrirang Abhyankar   Tao            tao=(Tao)ctx;
562*aad13602SShrirang Abhyankar   TAO_PDIPM      *pdipm = (TAO_PDIPM*)tao->data;
563*aad13602SShrirang Abhyankar   SNES           snes;
564*aad13602SShrirang Abhyankar   Vec            X,F,Y,W,G;
565*aad13602SShrirang Abhyankar   PetscInt       i,iter;
566*aad13602SShrirang Abhyankar   PetscReal      alpha_p=1.0,alpha_d=1.0,alpha[4];
567*aad13602SShrirang Abhyankar   PetscScalar    *Xarr,*z,*lambdai,dot;
568*aad13602SShrirang Abhyankar   const PetscScalar *dXarr,*dz,*dlambdai;
569*aad13602SShrirang Abhyankar   PetscScalar    *taosolarr;
570*aad13602SShrirang Abhyankar 
571*aad13602SShrirang Abhyankar   PetscFunctionBegin;
572*aad13602SShrirang Abhyankar   ierr = SNESLineSearchGetSNES(linesearch,&snes);CHKERRQ(ierr);
573*aad13602SShrirang Abhyankar   ierr = SNESGetIterationNumber(snes,&iter);CHKERRQ(ierr);
574*aad13602SShrirang Abhyankar 
575*aad13602SShrirang Abhyankar   ierr = SNESLineSearchSetReason(linesearch,SNES_LINESEARCH_SUCCEEDED);CHKERRQ(ierr);
576*aad13602SShrirang Abhyankar   ierr = SNESLineSearchGetVecs(linesearch,&X,&F,&Y,&W,&G);CHKERRQ(ierr);
577*aad13602SShrirang Abhyankar 
578*aad13602SShrirang Abhyankar   ierr = VecGetArray(X,&Xarr);CHKERRQ(ierr);
579*aad13602SShrirang Abhyankar   ierr = VecGetArrayRead(Y,&dXarr);CHKERRQ(ierr);
580*aad13602SShrirang Abhyankar   z  = Xarr + pdipm->off_z;
581*aad13602SShrirang Abhyankar   dz = dXarr + pdipm->off_z;
582*aad13602SShrirang Abhyankar   for (i=0; i < pdipm->nci; i++) {
583*aad13602SShrirang Abhyankar     if (z[i] - dz[i] < 0.0) {
584*aad13602SShrirang Abhyankar       alpha_p = PetscMin(alpha_p,0.9999*z[i]/dz[i]);
585*aad13602SShrirang Abhyankar     }
586*aad13602SShrirang Abhyankar   }
587*aad13602SShrirang Abhyankar 
588*aad13602SShrirang Abhyankar   lambdai  = Xarr + pdipm->off_lambdai;
589*aad13602SShrirang Abhyankar   dlambdai = dXarr + pdipm->off_lambdai;
590*aad13602SShrirang Abhyankar 
591*aad13602SShrirang Abhyankar   for (i=0; i<pdipm->nci; i++) {
592*aad13602SShrirang Abhyankar     if (lambdai[i] - dlambdai[i] < 0.0) {
593*aad13602SShrirang Abhyankar       alpha_d = PetscMin(0.9999*lambdai[i]/dlambdai[i],alpha_d);
594*aad13602SShrirang Abhyankar     }
595*aad13602SShrirang Abhyankar   }
596*aad13602SShrirang Abhyankar 
597*aad13602SShrirang Abhyankar   alpha[0] = alpha_p;
598*aad13602SShrirang Abhyankar   alpha[1] = alpha_d;
599*aad13602SShrirang Abhyankar   ierr = VecRestoreArrayRead(Y,&dXarr);CHKERRQ(ierr);
600*aad13602SShrirang Abhyankar   ierr = VecRestoreArray(X,&Xarr);CHKERRQ(ierr);
601*aad13602SShrirang Abhyankar 
602*aad13602SShrirang Abhyankar   /* alpha = min(alpha) over all processes */
603*aad13602SShrirang Abhyankar   ierr = MPI_Allreduce(alpha,alpha+2,2,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)tao));CHKERRQ(ierr);
604*aad13602SShrirang Abhyankar 
605*aad13602SShrirang Abhyankar   alpha_p = alpha[2];
606*aad13602SShrirang Abhyankar   alpha_d = alpha[3];
607*aad13602SShrirang Abhyankar 
608*aad13602SShrirang Abhyankar   ierr = VecGetArray(X,&Xarr);CHKERRQ(ierr);
609*aad13602SShrirang Abhyankar   ierr = VecGetArrayRead(Y,&dXarr);CHKERRQ(ierr);
610*aad13602SShrirang Abhyankar   for (i=0; i<pdipm->nx; i++) {
611*aad13602SShrirang Abhyankar     Xarr[i] = Xarr[i] - alpha_p * dXarr[i];
612*aad13602SShrirang Abhyankar   }
613*aad13602SShrirang Abhyankar 
614*aad13602SShrirang Abhyankar   for (i=0; i<pdipm->nce; i++) {
615*aad13602SShrirang Abhyankar     Xarr[i+pdipm->off_lambdae] = Xarr[i+pdipm->off_lambdae] - alpha_d * dXarr[i+pdipm->off_lambdae];
616*aad13602SShrirang Abhyankar   }
617*aad13602SShrirang Abhyankar 
618*aad13602SShrirang Abhyankar   for (i=0; i<pdipm->nci; i++) {
619*aad13602SShrirang Abhyankar     Xarr[i+pdipm->off_lambdai] = Xarr[i+pdipm->off_lambdai] - alpha_d * dXarr[i+pdipm->off_lambdai];
620*aad13602SShrirang Abhyankar   }
621*aad13602SShrirang Abhyankar 
622*aad13602SShrirang Abhyankar   for (i=0; i<pdipm->nci; i++) {
623*aad13602SShrirang Abhyankar     Xarr[i+pdipm->off_z] = Xarr[i+pdipm->off_z] - alpha_p * dXarr[i+pdipm->off_z];
624*aad13602SShrirang Abhyankar   }
625*aad13602SShrirang Abhyankar 
626*aad13602SShrirang Abhyankar   ierr = VecGetArray(tao->solution,&taosolarr);CHKERRQ(ierr);
627*aad13602SShrirang Abhyankar   ierr = PetscMemcpy(taosolarr,Xarr,pdipm->nx*sizeof(PetscScalar));CHKERRQ(ierr);
628*aad13602SShrirang Abhyankar   ierr = VecRestoreArray(tao->solution,&taosolarr);CHKERRQ(ierr);
629*aad13602SShrirang Abhyankar 
630*aad13602SShrirang Abhyankar 
631*aad13602SShrirang Abhyankar   ierr = VecRestoreArray(X,&Xarr);CHKERRQ(ierr);
632*aad13602SShrirang Abhyankar   ierr = VecRestoreArrayRead(Y,&dXarr);CHKERRQ(ierr);
633*aad13602SShrirang Abhyankar 
634*aad13602SShrirang Abhyankar   /* Evaluate F at X */
635*aad13602SShrirang Abhyankar   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
636*aad13602SShrirang Abhyankar   ierr = SNESLineSearchComputeNorms(linesearch);CHKERRQ(ierr); /* must call this func, do not know why */
637*aad13602SShrirang Abhyankar 
638*aad13602SShrirang Abhyankar   /* update mu = mu_update_factor * dot(z,lambdai)/pdipm->nci at updated X */
639*aad13602SShrirang Abhyankar   ierr = VecDot(pdipm->z,pdipm->lambdai,&dot);CHKERRQ(ierr);
640*aad13602SShrirang Abhyankar 
641*aad13602SShrirang Abhyankar   /* if (PetscAbsReal(pdipm->gradL) < 0.9*pdipm->mu)  */
642*aad13602SShrirang Abhyankar   pdipm->mu = pdipm->mu_update_factor * dot/pdipm->Nci;
643*aad13602SShrirang Abhyankar 
644*aad13602SShrirang Abhyankar   /* Update F; get tao->residual and tao->cnorm */
645*aad13602SShrirang Abhyankar   ierr = TaoSNESFunction_PDIPM(snes,X,F,(void*)tao);CHKERRQ(ierr);
646*aad13602SShrirang Abhyankar 
647*aad13602SShrirang Abhyankar   tao->niter++;
648*aad13602SShrirang Abhyankar   ierr = TaoLogConvergenceHistory(tao,pdipm->obj,tao->residual,tao->cnorm,tao->niter);CHKERRQ(ierr);
649*aad13602SShrirang Abhyankar   ierr = TaoMonitor(tao,tao->niter,pdipm->obj,tao->residual,tao->cnorm,pdipm->mu);CHKERRQ(ierr);
650*aad13602SShrirang Abhyankar 
651*aad13602SShrirang Abhyankar   ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
652*aad13602SShrirang Abhyankar   if (tao->reason) {
653*aad13602SShrirang Abhyankar     ierr = SNESSetConvergedReason(snes,SNES_CONVERGED_FNORM_ABS);CHKERRQ(ierr);
654*aad13602SShrirang Abhyankar   }
655*aad13602SShrirang Abhyankar   PetscFunctionReturn(0);
656*aad13602SShrirang Abhyankar }
657*aad13602SShrirang Abhyankar 
658*aad13602SShrirang Abhyankar /*
659*aad13602SShrirang Abhyankar    TaoSolve_PDIPM
660*aad13602SShrirang Abhyankar 
661*aad13602SShrirang Abhyankar    Input Parameter:
662*aad13602SShrirang Abhyankar    tao - TAO context
663*aad13602SShrirang Abhyankar 
664*aad13602SShrirang Abhyankar    Output Parameter:
665*aad13602SShrirang Abhyankar    tao - TAO context
666*aad13602SShrirang Abhyankar */
667*aad13602SShrirang Abhyankar PetscErrorCode TaoSolve_PDIPM(Tao tao)
668*aad13602SShrirang Abhyankar {
669*aad13602SShrirang Abhyankar   PetscErrorCode     ierr;
670*aad13602SShrirang Abhyankar   TAO_PDIPM          *pdipm = (TAO_PDIPM*)tao->data;
671*aad13602SShrirang Abhyankar   SNESLineSearch     linesearch;  /* SNESLineSearch context */
672*aad13602SShrirang Abhyankar   Vec                dummy;
673*aad13602SShrirang Abhyankar 
674*aad13602SShrirang Abhyankar   PetscFunctionBegin;
675*aad13602SShrirang Abhyankar   /* Initialize all variables */
676*aad13602SShrirang Abhyankar   ierr = TaoPDIPMInitializeSolution(tao);CHKERRQ(ierr);
677*aad13602SShrirang Abhyankar 
678*aad13602SShrirang Abhyankar   /* Set linesearch */
679*aad13602SShrirang Abhyankar   ierr = SNESGetLineSearch(pdipm->snes,&linesearch);CHKERRQ(ierr);
680*aad13602SShrirang Abhyankar   ierr = SNESLineSearchSetType(linesearch,SNESLINESEARCHSHELL);CHKERRQ(ierr);
681*aad13602SShrirang Abhyankar   ierr = SNESLineSearchShellSetUserFunc(linesearch,PDIPMLineSearch,tao);CHKERRQ(ierr);
682*aad13602SShrirang Abhyankar   ierr = SNESLineSearchSetFromOptions(linesearch);CHKERRQ(ierr);
683*aad13602SShrirang Abhyankar 
684*aad13602SShrirang Abhyankar   tao->reason = TAO_CONTINUE_ITERATING;
685*aad13602SShrirang Abhyankar 
686*aad13602SShrirang Abhyankar   /* -tao_monitor for iteration 0 and check convergence */
687*aad13602SShrirang Abhyankar   ierr = VecDuplicate(pdipm->X,&dummy);CHKERRQ(ierr);
688*aad13602SShrirang Abhyankar   ierr = TaoSNESFunction_PDIPM(pdipm->snes,pdipm->X,dummy,(void*)tao);CHKERRQ(ierr);
689*aad13602SShrirang Abhyankar 
690*aad13602SShrirang Abhyankar   ierr = TaoLogConvergenceHistory(tao,pdipm->obj,tao->residual,tao->cnorm,tao->niter);CHKERRQ(ierr);
691*aad13602SShrirang Abhyankar   ierr = TaoMonitor(tao,tao->niter,pdipm->obj,tao->residual,tao->cnorm,pdipm->mu);CHKERRQ(ierr);
692*aad13602SShrirang Abhyankar   ierr = VecDestroy(&dummy);CHKERRQ(ierr);
693*aad13602SShrirang Abhyankar   ierr = (*tao->ops->convergencetest)(tao,tao->cnvP);CHKERRQ(ierr);
694*aad13602SShrirang Abhyankar   if (tao->reason) {
695*aad13602SShrirang Abhyankar     ierr = SNESSetConvergedReason(pdipm->snes,SNES_CONVERGED_FNORM_ABS);CHKERRQ(ierr);
696*aad13602SShrirang Abhyankar   }
697*aad13602SShrirang Abhyankar 
698*aad13602SShrirang Abhyankar   while (tao->reason == TAO_CONTINUE_ITERATING) {
699*aad13602SShrirang Abhyankar     SNESConvergedReason reason;
700*aad13602SShrirang Abhyankar     ierr = SNESSolve(pdipm->snes,NULL,pdipm->X);CHKERRQ(ierr);
701*aad13602SShrirang Abhyankar 
702*aad13602SShrirang Abhyankar     /* Check SNES convergence */
703*aad13602SShrirang Abhyankar     ierr = SNESGetConvergedReason(pdipm->snes,&reason);CHKERRQ(ierr);
704*aad13602SShrirang Abhyankar     if (reason < 0) {
705*aad13602SShrirang Abhyankar       ierr = PetscPrintf(PETSC_COMM_WORLD,"SNES solve did not converged due to reason %D\n",reason);CHKERRQ(ierr);
706*aad13602SShrirang Abhyankar     }
707*aad13602SShrirang Abhyankar 
708*aad13602SShrirang Abhyankar     /* Check TAO convergence */
709*aad13602SShrirang Abhyankar     if (PetscIsInfOrNanReal(pdipm->obj)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"User-provided compute function generated Inf or NaN");
710*aad13602SShrirang Abhyankar   }
711*aad13602SShrirang Abhyankar   PetscFunctionReturn(0);
712*aad13602SShrirang Abhyankar }
713*aad13602SShrirang Abhyankar 
714*aad13602SShrirang Abhyankar /*
715*aad13602SShrirang Abhyankar    TaoSetup_PDIPM - Sets up tao and pdipm
716*aad13602SShrirang Abhyankar 
717*aad13602SShrirang Abhyankar    Input Parameter:
718*aad13602SShrirang Abhyankar    tao - TAO object
719*aad13602SShrirang Abhyankar 
720*aad13602SShrirang Abhyankar    Output:   pdipm - initialized object
721*aad13602SShrirang Abhyankar */
722*aad13602SShrirang Abhyankar PetscErrorCode TaoSetup_PDIPM(Tao tao)
723*aad13602SShrirang Abhyankar {
724*aad13602SShrirang Abhyankar   TAO_PDIPM      *pdipm = (TAO_PDIPM*)tao->data;
725*aad13602SShrirang Abhyankar   PetscErrorCode ierr;
726*aad13602SShrirang Abhyankar   MPI_Comm       comm;
727*aad13602SShrirang Abhyankar   PetscMPIInt    rank,size;
728*aad13602SShrirang Abhyankar   PetscInt       row,col,Jcrstart,Jcrend,k,tmp,nc,proc,*nh_all,*ng_all;
729*aad13602SShrirang Abhyankar   PetscInt       offset,*xa,*xb,i,j,rstart,rend;
730*aad13602SShrirang Abhyankar   PetscScalar    one=1.0,neg_one=-1.0,*Xarr;
731*aad13602SShrirang Abhyankar   const PetscInt    *cols,*rranges,*cranges,*aj,*ranges;
732*aad13602SShrirang Abhyankar   const PetscScalar *aa;
733*aad13602SShrirang Abhyankar   Mat            J,jac_equality_trans,jac_inequality_trans;
734*aad13602SShrirang Abhyankar   Mat            Jce_xfixed_trans,Jci_xb_trans;
735*aad13602SShrirang Abhyankar   PetscInt       *dnz,*onz,rjstart,nx_all,*nce_all,*Jranges,cols1[2];
736*aad13602SShrirang Abhyankar 
737*aad13602SShrirang Abhyankar   PetscFunctionBegin;
738*aad13602SShrirang Abhyankar   ierr = PetscObjectGetComm((PetscObject)tao,&comm);CHKERRQ(ierr);
739*aad13602SShrirang Abhyankar   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
740*aad13602SShrirang Abhyankar   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
741*aad13602SShrirang Abhyankar 
742*aad13602SShrirang Abhyankar   /* (1) Setup Bounds and create Tao vectors */
743*aad13602SShrirang Abhyankar   ierr = TaoPDIPMSetUpBounds(tao);CHKERRQ(ierr);
744*aad13602SShrirang Abhyankar 
745*aad13602SShrirang Abhyankar   if (!tao->gradient) {
746*aad13602SShrirang Abhyankar     ierr = VecDuplicate(tao->solution,&tao->gradient);CHKERRQ(ierr);
747*aad13602SShrirang Abhyankar     ierr = VecDuplicate(tao->solution,&tao->stepdirection);CHKERRQ(ierr);
748*aad13602SShrirang Abhyankar   }
749*aad13602SShrirang Abhyankar 
750*aad13602SShrirang Abhyankar   /* (2) Get sizes */
751*aad13602SShrirang Abhyankar   /* Size of vector x - This is set by TaoSetInitia√lVector */
752*aad13602SShrirang Abhyankar   ierr = VecGetSize(tao->solution,&pdipm->Nx);CHKERRQ(ierr);
753*aad13602SShrirang Abhyankar   ierr = VecGetLocalSize(tao->solution,&pdipm->nx);CHKERRQ(ierr);
754*aad13602SShrirang Abhyankar 
755*aad13602SShrirang Abhyankar   /* Size of equality constraints and vectors */
756*aad13602SShrirang Abhyankar   if (tao->constraints_equality) {
757*aad13602SShrirang Abhyankar     ierr = VecGetSize(tao->constraints_equality,&pdipm->Ng);CHKERRQ(ierr);
758*aad13602SShrirang Abhyankar     ierr = VecGetLocalSize(tao->constraints_equality,&pdipm->ng);CHKERRQ(ierr);
759*aad13602SShrirang Abhyankar   } else {
760*aad13602SShrirang Abhyankar     pdipm->ng = pdipm->Ng = 0;
761*aad13602SShrirang Abhyankar   }
762*aad13602SShrirang Abhyankar 
763*aad13602SShrirang Abhyankar   pdipm->nce = pdipm->ng + pdipm->nxfixed;
764*aad13602SShrirang Abhyankar   pdipm->Nce = pdipm->Ng + pdipm->Nxfixed;
765*aad13602SShrirang Abhyankar 
766*aad13602SShrirang Abhyankar   /* Size of inequality constraints and vectors */
767*aad13602SShrirang Abhyankar   if (tao->constraints_inequality) {
768*aad13602SShrirang Abhyankar     ierr = VecGetSize(tao->constraints_inequality,&pdipm->Nh);CHKERRQ(ierr);
769*aad13602SShrirang Abhyankar     ierr = VecGetLocalSize(tao->constraints_inequality,&pdipm->nh);CHKERRQ(ierr);
770*aad13602SShrirang Abhyankar   } else {
771*aad13602SShrirang Abhyankar     pdipm->nh = pdipm->Nh = 0;
772*aad13602SShrirang Abhyankar   }
773*aad13602SShrirang Abhyankar 
774*aad13602SShrirang Abhyankar   pdipm->nci = pdipm->nh + pdipm->nxlb + pdipm->nxub + 2*pdipm->nxbox;
775*aad13602SShrirang Abhyankar   pdipm->Nci = pdipm->Nh + pdipm->Nxlb + pdipm->Nxub + 2*pdipm->Nxbox;
776*aad13602SShrirang Abhyankar 
777*aad13602SShrirang Abhyankar   /* Full size of the KKT system to be solved */
778*aad13602SShrirang Abhyankar   pdipm->n = pdipm->nx + pdipm->nce + 2*pdipm->nci;
779*aad13602SShrirang Abhyankar   pdipm->N = pdipm->Nx + pdipm->Nce + 2*pdipm->Nci;
780*aad13602SShrirang Abhyankar 
781*aad13602SShrirang Abhyankar   /* list below to TaoView_PDIPM()? */
782*aad13602SShrirang Abhyankar   /* ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] nce %d = ng %d + nxfixed %d\n",rank,pdipm->nce,pdipm->ng,pdipm->nxfixed); */
783*aad13602SShrirang Abhyankar   /* ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] nci %d = nh %d + nxlb %d + nxub %d + 2*nxbox %d\n",rank,pdipm->nci,pdipm->nh,pdipm->nxlb,pdipm->nxub,pdipm->nxbox); */
784*aad13602SShrirang Abhyankar   /* ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] n %d = nx %d + nce %d + 2*nci %d\n",rank,pdipm->n,pdipm->nx,pdipm->nce,pdipm->nci); */
785*aad13602SShrirang Abhyankar 
786*aad13602SShrirang Abhyankar   /* (3) Offsets for subvectors */
787*aad13602SShrirang Abhyankar   pdipm->off_lambdae = pdipm->nx;
788*aad13602SShrirang Abhyankar   pdipm->off_lambdai = pdipm->off_lambdae + pdipm->nce;
789*aad13602SShrirang Abhyankar   pdipm->off_z       = pdipm->off_lambdai + pdipm->nci;
790*aad13602SShrirang Abhyankar 
791*aad13602SShrirang Abhyankar   /* (4) Create vectors and subvectors */
792*aad13602SShrirang Abhyankar   /* Ce and Ci vectors */
793*aad13602SShrirang Abhyankar   ierr = VecCreate(comm,&pdipm->ce);CHKERRQ(ierr);
794*aad13602SShrirang Abhyankar   ierr = VecSetSizes(pdipm->ce,pdipm->nce,pdipm->Nce);CHKERRQ(ierr);
795*aad13602SShrirang Abhyankar   ierr = VecSetFromOptions(pdipm->ce);CHKERRQ(ierr);
796*aad13602SShrirang Abhyankar 
797*aad13602SShrirang Abhyankar   ierr = VecCreate(comm,&pdipm->ci);CHKERRQ(ierr);
798*aad13602SShrirang Abhyankar   ierr = VecSetSizes(pdipm->ci,pdipm->nci,pdipm->Nci);CHKERRQ(ierr);
799*aad13602SShrirang Abhyankar   ierr = VecSetFromOptions(pdipm->ci);CHKERRQ(ierr);
800*aad13602SShrirang Abhyankar 
801*aad13602SShrirang Abhyankar   /* X=[x; lambdae; lambdai; z] for the big KKT system */
802*aad13602SShrirang Abhyankar   ierr = VecCreate(comm,&pdipm->X);CHKERRQ(ierr);
803*aad13602SShrirang Abhyankar   ierr = VecSetSizes(pdipm->X,pdipm->n,pdipm->N);CHKERRQ(ierr);
804*aad13602SShrirang Abhyankar   ierr = VecSetFromOptions(pdipm->X);CHKERRQ(ierr);
805*aad13602SShrirang Abhyankar 
806*aad13602SShrirang Abhyankar   /* Subvectors; they share local arrays with X */
807*aad13602SShrirang Abhyankar   ierr = VecGetArray(pdipm->X,&Xarr);CHKERRQ(ierr);
808*aad13602SShrirang Abhyankar   /* x shares local array with X.x */
809*aad13602SShrirang Abhyankar   if (pdipm->Nx) {
810*aad13602SShrirang Abhyankar     ierr = VecCreateMPIWithArray(comm,1,pdipm->nx,pdipm->Nx,Xarr,&pdipm->x);CHKERRQ(ierr);
811*aad13602SShrirang Abhyankar   }
812*aad13602SShrirang Abhyankar 
813*aad13602SShrirang Abhyankar   /* lambdae shares local array with X.lambdae */
814*aad13602SShrirang Abhyankar   if (pdipm->Nce) {
815*aad13602SShrirang Abhyankar     ierr = VecCreateMPIWithArray(comm,1,pdipm->nce,pdipm->Nce,Xarr+pdipm->off_lambdae,&pdipm->lambdae);CHKERRQ(ierr);
816*aad13602SShrirang Abhyankar   }
817*aad13602SShrirang Abhyankar 
818*aad13602SShrirang Abhyankar   /* tao->DE shares local array with X.lambdae_g */
819*aad13602SShrirang Abhyankar   if (pdipm->Ng) {
820*aad13602SShrirang Abhyankar     ierr = VecCreateMPIWithArray(comm,1,pdipm->ng,pdipm->Ng,Xarr+pdipm->off_lambdae,&tao->DE);CHKERRQ(ierr);
821*aad13602SShrirang Abhyankar 
822*aad13602SShrirang Abhyankar     ierr = VecCreate(comm,&pdipm->lambdae_xfixed);CHKERRQ(ierr);
823*aad13602SShrirang Abhyankar     ierr = VecSetSizes(pdipm->lambdae_xfixed,pdipm->nxfixed,PETSC_DECIDE);CHKERRQ(ierr);
824*aad13602SShrirang Abhyankar     ierr = VecSetFromOptions(pdipm->lambdae_xfixed);CHKERRQ(ierr);
825*aad13602SShrirang Abhyankar   }
826*aad13602SShrirang Abhyankar 
827*aad13602SShrirang Abhyankar   if (pdipm->Nci) {
828*aad13602SShrirang Abhyankar     /* lambdai shares local array with X.lambdai */
829*aad13602SShrirang Abhyankar     ierr = VecCreateMPIWithArray(comm,1,pdipm->nci,pdipm->Nci,Xarr+pdipm->off_lambdai,&pdipm->lambdai);CHKERRQ(ierr);
830*aad13602SShrirang Abhyankar 
831*aad13602SShrirang Abhyankar     /* z for slack variables; it shares local array with X.z */
832*aad13602SShrirang Abhyankar     ierr = VecCreateMPIWithArray(comm,1,pdipm->nci,pdipm->Nci,Xarr+pdipm->off_z,&pdipm->z);CHKERRQ(ierr);
833*aad13602SShrirang Abhyankar   }
834*aad13602SShrirang Abhyankar 
835*aad13602SShrirang Abhyankar   /* tao->DI which shares local array with X.lambdai_h */
836*aad13602SShrirang Abhyankar   if (pdipm->Nh) {
837*aad13602SShrirang Abhyankar     ierr = VecCreateMPIWithArray(comm,1,pdipm->nh,pdipm->Nh,Xarr+pdipm->off_lambdai,&tao->DI);CHKERRQ(ierr);
838*aad13602SShrirang Abhyankar   }
839*aad13602SShrirang Abhyankar 
840*aad13602SShrirang Abhyankar   ierr = VecCreate(comm,&pdipm->lambdai_xb);CHKERRQ(ierr);
841*aad13602SShrirang Abhyankar   ierr = VecSetSizes(pdipm->lambdai_xb,(pdipm->nci - pdipm->nh),PETSC_DECIDE);CHKERRQ(ierr);
842*aad13602SShrirang Abhyankar   ierr = VecSetFromOptions(pdipm->lambdai_xb);CHKERRQ(ierr);
843*aad13602SShrirang Abhyankar 
844*aad13602SShrirang Abhyankar   ierr = VecRestoreArray(pdipm->X,&Xarr);CHKERRQ(ierr);
845*aad13602SShrirang Abhyankar 
846*aad13602SShrirang Abhyankar   /* (5) Create Jacobians Jce_xfixed and Jci */
847*aad13602SShrirang Abhyankar   /* (5.1) PDIPM Jacobian of equality bounds cebound(x) = J_nxfixed */
848*aad13602SShrirang Abhyankar   if (pdipm->Nxfixed) {
849*aad13602SShrirang Abhyankar     /* Create Jce_xfixed */
850*aad13602SShrirang Abhyankar     ierr = MatCreate(comm,&pdipm->Jce_xfixed);CHKERRQ(ierr);
851*aad13602SShrirang Abhyankar     ierr = MatSetSizes(pdipm->Jce_xfixed,pdipm->nxfixed,pdipm->nx,PETSC_DECIDE,pdipm->Nx);CHKERRQ(ierr);
852*aad13602SShrirang Abhyankar     ierr = MatSetFromOptions(pdipm->Jce_xfixed);CHKERRQ(ierr);
853*aad13602SShrirang Abhyankar     ierr = MatSeqAIJSetPreallocation(pdipm->Jce_xfixed,1,NULL);CHKERRQ(ierr);
854*aad13602SShrirang Abhyankar     ierr = MatMPIAIJSetPreallocation(pdipm->Jce_xfixed,1,NULL,1,NULL);CHKERRQ(ierr);
855*aad13602SShrirang Abhyankar 
856*aad13602SShrirang Abhyankar     ierr = MatGetOwnershipRange(pdipm->Jce_xfixed,&Jcrstart,&Jcrend);CHKERRQ(ierr);
857*aad13602SShrirang Abhyankar     ierr = ISGetIndices(pdipm->isxfixed,&cols);CHKERRQ(ierr);
858*aad13602SShrirang Abhyankar     k = 0;
859*aad13602SShrirang Abhyankar     for (row = Jcrstart; row < Jcrend; row++) {
860*aad13602SShrirang Abhyankar       ierr = MatSetValues(pdipm->Jce_xfixed,1,&row,1,cols+k,&one,INSERT_VALUES);CHKERRQ(ierr);
861*aad13602SShrirang Abhyankar       k++;
862*aad13602SShrirang Abhyankar     }
863*aad13602SShrirang Abhyankar     ierr = ISRestoreIndices(pdipm->isxfixed, &cols);CHKERRQ(ierr);
864*aad13602SShrirang Abhyankar     ierr = MatAssemblyBegin(pdipm->Jce_xfixed,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
865*aad13602SShrirang Abhyankar     ierr = MatAssemblyEnd(pdipm->Jce_xfixed,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
866*aad13602SShrirang Abhyankar   }
867*aad13602SShrirang Abhyankar 
868*aad13602SShrirang Abhyankar   /* (5.2) PDIPM inequality Jacobian Jci = [tao->jacobian_inequality; ...] */
869*aad13602SShrirang Abhyankar   ierr = MatCreate(comm,&pdipm->Jci_xb);CHKERRQ(ierr);
870*aad13602SShrirang Abhyankar   ierr = MatSetSizes(pdipm->Jci_xb,pdipm->nci-pdipm->nh,pdipm->nx,PETSC_DECIDE,pdipm->Nx);CHKERRQ(ierr);
871*aad13602SShrirang Abhyankar   ierr = MatSetFromOptions(pdipm->Jci_xb);CHKERRQ(ierr);
872*aad13602SShrirang Abhyankar   ierr = MatSeqAIJSetPreallocation(pdipm->Jci_xb,1,NULL);CHKERRQ(ierr);
873*aad13602SShrirang Abhyankar   ierr = MatMPIAIJSetPreallocation(pdipm->Jci_xb,1,NULL,1,NULL);CHKERRQ(ierr);
874*aad13602SShrirang Abhyankar 
875*aad13602SShrirang Abhyankar   ierr = MatGetOwnershipRange(pdipm->Jci_xb,&Jcrstart,&Jcrend);CHKERRQ(ierr);
876*aad13602SShrirang Abhyankar   offset = Jcrstart;
877*aad13602SShrirang Abhyankar   if (pdipm->Nxub) {
878*aad13602SShrirang Abhyankar     /* Add xub to Jci_xb */
879*aad13602SShrirang Abhyankar     ierr = ISGetIndices(pdipm->isxub,&cols);CHKERRQ(ierr);
880*aad13602SShrirang Abhyankar     k = 0;
881*aad13602SShrirang Abhyankar     for (row = offset; row < offset + pdipm->nxub; row++) {
882*aad13602SShrirang Abhyankar       ierr = MatSetValues(pdipm->Jci_xb,1,&row,1,cols+k,&neg_one,INSERT_VALUES);CHKERRQ(ierr);
883*aad13602SShrirang Abhyankar       k++;
884*aad13602SShrirang Abhyankar     }
885*aad13602SShrirang Abhyankar     ierr = ISRestoreIndices(pdipm->isxub, &cols);CHKERRQ(ierr);
886*aad13602SShrirang Abhyankar   }
887*aad13602SShrirang Abhyankar 
888*aad13602SShrirang Abhyankar   if (pdipm->Nxlb) {
889*aad13602SShrirang Abhyankar     /* Add xlb to Jci_xb */
890*aad13602SShrirang Abhyankar     ierr = ISGetIndices(pdipm->isxlb,&cols);CHKERRQ(ierr);
891*aad13602SShrirang Abhyankar     k = 0;
892*aad13602SShrirang Abhyankar     offset += pdipm->nxub;
893*aad13602SShrirang Abhyankar     for (row = offset; row < offset + pdipm->nxlb; row++) {
894*aad13602SShrirang Abhyankar       ierr = MatSetValues(pdipm->Jci_xb,1,&row,1,cols+k,&one,INSERT_VALUES);CHKERRQ(ierr);
895*aad13602SShrirang Abhyankar       k++;
896*aad13602SShrirang Abhyankar     }
897*aad13602SShrirang Abhyankar     ierr = ISRestoreIndices(pdipm->isxlb, &cols);CHKERRQ(ierr);
898*aad13602SShrirang Abhyankar   }
899*aad13602SShrirang Abhyankar 
900*aad13602SShrirang Abhyankar   /* Add xbox to Jci_xb */
901*aad13602SShrirang Abhyankar   if (pdipm->Nxbox) {
902*aad13602SShrirang Abhyankar     ierr = ISGetIndices(pdipm->isxbox,&cols);CHKERRQ(ierr);
903*aad13602SShrirang Abhyankar     k = 0;
904*aad13602SShrirang Abhyankar     offset += pdipm->nxlb;
905*aad13602SShrirang Abhyankar     for (row = offset; row < offset + pdipm->nxbox; row++) {
906*aad13602SShrirang Abhyankar       ierr = MatSetValues(pdipm->Jci_xb,1,&row,1,cols+k,&neg_one,INSERT_VALUES);CHKERRQ(ierr);
907*aad13602SShrirang Abhyankar       tmp = row + pdipm->nxbox;
908*aad13602SShrirang Abhyankar       ierr = MatSetValues(pdipm->Jci_xb,1,&tmp,1,cols+k,&one,INSERT_VALUES);CHKERRQ(ierr);
909*aad13602SShrirang Abhyankar       k++;
910*aad13602SShrirang Abhyankar     }
911*aad13602SShrirang Abhyankar     ierr = ISRestoreIndices(pdipm->isxbox, &cols);CHKERRQ(ierr);
912*aad13602SShrirang Abhyankar   }
913*aad13602SShrirang Abhyankar 
914*aad13602SShrirang Abhyankar   ierr = MatAssemblyBegin(pdipm->Jci_xb,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
915*aad13602SShrirang Abhyankar   ierr = MatAssemblyEnd(pdipm->Jci_xb,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
916*aad13602SShrirang Abhyankar   /* ierr = MatView(pdipm->Jci_xb,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */
917*aad13602SShrirang Abhyankar 
918*aad13602SShrirang Abhyankar   /* (6) Set up ISs for PC Fieldsplit */
919*aad13602SShrirang Abhyankar   if (pdipm->solve_reduced_kkt) {
920*aad13602SShrirang Abhyankar     ierr = PetscMalloc2(pdipm->nx+pdipm->nce,&xa,2*pdipm->nci,&xb);CHKERRQ(ierr);
921*aad13602SShrirang Abhyankar     for(i=0; i < pdipm->nx + pdipm->nce; i++) xa[i] = i;
922*aad13602SShrirang Abhyankar     for(i=0; i < 2*pdipm->nci; i++) xb[i] = pdipm->off_lambdai + i;
923*aad13602SShrirang Abhyankar 
924*aad13602SShrirang Abhyankar     ierr = ISCreateGeneral(comm,pdipm->nx+pdipm->nce,xa,PETSC_OWN_POINTER,&pdipm->is1);CHKERRQ(ierr);
925*aad13602SShrirang Abhyankar     ierr = ISCreateGeneral(comm,2*pdipm->nci,xb,PETSC_OWN_POINTER,&pdipm->is2);CHKERRQ(ierr);
926*aad13602SShrirang Abhyankar   }
927*aad13602SShrirang Abhyankar 
928*aad13602SShrirang Abhyankar   /* (7) Gather offsets from all processes */
929*aad13602SShrirang Abhyankar   ierr = PetscMalloc1(size,&pdipm->nce_all);CHKERRQ(ierr);
930*aad13602SShrirang Abhyankar 
931*aad13602SShrirang Abhyankar   /* Get rstart of KKT matrix */
932*aad13602SShrirang Abhyankar   ierr = MPI_Scan(&pdipm->n,&rstart,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr);
933*aad13602SShrirang Abhyankar   rstart -= pdipm->n;
934*aad13602SShrirang Abhyankar 
935*aad13602SShrirang Abhyankar   ierr = MPI_Allgather(&pdipm->nce,1,MPIU_INT,pdipm->nce_all,1,MPIU_INT,comm);CHKERRQ(ierr);
936*aad13602SShrirang Abhyankar 
937*aad13602SShrirang Abhyankar   ierr = PetscMalloc3(size,&ng_all,size,&nh_all,size,&Jranges);CHKERRQ(ierr);
938*aad13602SShrirang Abhyankar   ierr = MPI_Allgather(&rstart,1,MPIU_INT,Jranges,1,MPIU_INT,comm);CHKERRQ(ierr);
939*aad13602SShrirang Abhyankar   ierr = MPI_Allgather(&pdipm->nh,1,MPIU_INT,nh_all,1,MPIU_INT,comm);CHKERRQ(ierr);
940*aad13602SShrirang Abhyankar   ierr = MPI_Allgather(&pdipm->ng,1,MPIU_INT,ng_all,1,MPIU_INT,comm);CHKERRQ(ierr);
941*aad13602SShrirang Abhyankar 
942*aad13602SShrirang Abhyankar   ierr = MatGetOwnershipRanges(tao->hessian,&rranges);CHKERRQ(ierr);
943*aad13602SShrirang Abhyankar   ierr = MatGetOwnershipRangesColumn(tao->hessian,&cranges);CHKERRQ(ierr);
944*aad13602SShrirang Abhyankar 
945*aad13602SShrirang Abhyankar   if (pdipm->Ng) {
946*aad13602SShrirang Abhyankar     ierr = TaoComputeJacobianEquality(tao,tao->solution,tao->jacobian_equality,tao->jacobian_equality_pre);CHKERRQ(ierr);
947*aad13602SShrirang Abhyankar     ierr = MatTranspose(tao->jacobian_equality,MAT_INITIAL_MATRIX,&pdipm->jac_equality_trans);CHKERRQ(ierr);
948*aad13602SShrirang Abhyankar   }
949*aad13602SShrirang Abhyankar   if (pdipm->Nh) {
950*aad13602SShrirang Abhyankar     ierr = TaoComputeJacobianInequality(tao,tao->solution,tao->jacobian_inequality,tao->jacobian_inequality_pre);CHKERRQ(ierr);
951*aad13602SShrirang Abhyankar     ierr = MatTranspose(tao->jacobian_inequality,MAT_INITIAL_MATRIX,&pdipm->jac_inequality_trans);CHKERRQ(ierr);
952*aad13602SShrirang Abhyankar   }
953*aad13602SShrirang Abhyankar 
954*aad13602SShrirang Abhyankar   /* Count dnz,onz for preallocation of KKT matrix */
955*aad13602SShrirang Abhyankar   jac_equality_trans   = pdipm->jac_equality_trans;
956*aad13602SShrirang Abhyankar   jac_inequality_trans = pdipm->jac_inequality_trans;
957*aad13602SShrirang Abhyankar   nce_all = pdipm->nce_all;
958*aad13602SShrirang Abhyankar 
959*aad13602SShrirang Abhyankar   if (pdipm->Nxfixed) {
960*aad13602SShrirang Abhyankar     ierr = MatTranspose(pdipm->Jce_xfixed,MAT_INITIAL_MATRIX,&Jce_xfixed_trans);CHKERRQ(ierr);
961*aad13602SShrirang Abhyankar   }
962*aad13602SShrirang Abhyankar   ierr = MatTranspose(pdipm->Jci_xb,MAT_INITIAL_MATRIX,&Jci_xb_trans);CHKERRQ(ierr);
963*aad13602SShrirang Abhyankar 
964*aad13602SShrirang Abhyankar   ierr = MatPreallocateInitialize(comm,pdipm->n,pdipm->n,dnz,onz);CHKERRQ(ierr);
965*aad13602SShrirang Abhyankar 
966*aad13602SShrirang Abhyankar   /* 1st row block of KKT matrix: [Wxx; gradCe'; -gradCi'; 0] */
967*aad13602SShrirang Abhyankar   ierr = TaoPDIPMEvaluateFunctionsAndJacobians(tao,pdipm->x);CHKERRQ(ierr);
968*aad13602SShrirang Abhyankar   ierr = TaoComputeHessian(tao,tao->solution,tao->hessian,tao->hessian_pre);CHKERRQ(ierr);
969*aad13602SShrirang Abhyankar 
970*aad13602SShrirang Abhyankar   /* Insert tao->hessian */
971*aad13602SShrirang Abhyankar   ierr = MatGetOwnershipRange(tao->hessian,&rjstart,NULL);CHKERRQ(ierr);
972*aad13602SShrirang Abhyankar   for (i=0; i<pdipm->nx; i++){
973*aad13602SShrirang Abhyankar     row = rstart + i;
974*aad13602SShrirang Abhyankar 
975*aad13602SShrirang Abhyankar     ierr = MatGetRow(tao->hessian,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr);
976*aad13602SShrirang Abhyankar     proc = 0;
977*aad13602SShrirang Abhyankar     for (j=0; j < nc; j++) {
978*aad13602SShrirang Abhyankar       while (aj[j] >= cranges[proc+1]) proc++;
979*aad13602SShrirang Abhyankar       col = aj[j] - cranges[proc] + Jranges[proc];
980*aad13602SShrirang Abhyankar       ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr);
981*aad13602SShrirang Abhyankar     }
982*aad13602SShrirang Abhyankar     ierr = MatRestoreRow(tao->hessian,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr);
983*aad13602SShrirang Abhyankar 
984*aad13602SShrirang Abhyankar     if(pdipm->ng) {
985*aad13602SShrirang Abhyankar       /* Insert grad g' */
986*aad13602SShrirang Abhyankar       ierr = MatGetRow(jac_equality_trans,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr);
987*aad13602SShrirang Abhyankar       ierr = MatGetOwnershipRanges(tao->jacobian_equality,&ranges);CHKERRQ(ierr);
988*aad13602SShrirang Abhyankar       proc = 0;
989*aad13602SShrirang Abhyankar       for (j=0; j < nc; j++) {
990*aad13602SShrirang Abhyankar         /* find row ownership of */
991*aad13602SShrirang Abhyankar         while (aj[j] >= ranges[proc+1]) proc++;
992*aad13602SShrirang Abhyankar         nx_all = rranges[proc+1] - rranges[proc];
993*aad13602SShrirang Abhyankar         col = aj[j] - ranges[proc] + Jranges[proc] + nx_all;
994*aad13602SShrirang Abhyankar         ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr);
995*aad13602SShrirang Abhyankar       }
996*aad13602SShrirang Abhyankar       ierr = MatRestoreRow(jac_equality_trans,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr);
997*aad13602SShrirang Abhyankar     }
998*aad13602SShrirang Abhyankar 
999*aad13602SShrirang Abhyankar     /* Insert Jce_xfixed^T' */
1000*aad13602SShrirang Abhyankar     if (pdipm->nxfixed) {
1001*aad13602SShrirang Abhyankar       ierr = MatGetRow(Jce_xfixed_trans,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr);
1002*aad13602SShrirang Abhyankar       ierr = MatGetOwnershipRanges(pdipm->Jce_xfixed,&ranges);CHKERRQ(ierr);
1003*aad13602SShrirang Abhyankar       proc = 0;
1004*aad13602SShrirang Abhyankar       for (j=0; j < nc; j++) {
1005*aad13602SShrirang Abhyankar         /* find row ownership of */
1006*aad13602SShrirang Abhyankar         while (aj[j] >= ranges[proc+1]) proc++;
1007*aad13602SShrirang Abhyankar         nx_all = rranges[proc+1] - rranges[proc];
1008*aad13602SShrirang Abhyankar         col = aj[j] - ranges[proc] + Jranges[proc] + nx_all + ng_all[proc];
1009*aad13602SShrirang Abhyankar         ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr);
1010*aad13602SShrirang Abhyankar       }
1011*aad13602SShrirang Abhyankar       ierr = MatRestoreRow(Jce_xfixed_trans,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr);
1012*aad13602SShrirang Abhyankar     }
1013*aad13602SShrirang Abhyankar 
1014*aad13602SShrirang Abhyankar     if(pdipm->nh) {
1015*aad13602SShrirang Abhyankar       /* Insert -grad h' */
1016*aad13602SShrirang Abhyankar       ierr = MatGetRow(jac_inequality_trans,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr);
1017*aad13602SShrirang Abhyankar       ierr = MatGetOwnershipRanges(tao->jacobian_inequality,&ranges);CHKERRQ(ierr);
1018*aad13602SShrirang Abhyankar       proc = 0;
1019*aad13602SShrirang Abhyankar       for (j=0; j < nc; j++) {
1020*aad13602SShrirang Abhyankar         /* find row ownership of */
1021*aad13602SShrirang Abhyankar         while (aj[j] >= ranges[proc+1]) proc++;
1022*aad13602SShrirang Abhyankar         nx_all = rranges[proc+1] - rranges[proc];
1023*aad13602SShrirang Abhyankar         col = aj[j] - ranges[proc] + Jranges[proc] + nx_all + nce_all[proc];
1024*aad13602SShrirang Abhyankar         ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr);
1025*aad13602SShrirang Abhyankar       }
1026*aad13602SShrirang Abhyankar       ierr = MatRestoreRow(jac_inequality_trans,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr);
1027*aad13602SShrirang Abhyankar     }
1028*aad13602SShrirang Abhyankar 
1029*aad13602SShrirang Abhyankar     /* Insert Jci_xb^T' */
1030*aad13602SShrirang Abhyankar     ierr = MatGetRow(Jci_xb_trans,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr);
1031*aad13602SShrirang Abhyankar     ierr = MatGetOwnershipRanges(pdipm->Jci_xb,&ranges);CHKERRQ(ierr);
1032*aad13602SShrirang Abhyankar     proc = 0;
1033*aad13602SShrirang Abhyankar     for (j=0; j < nc; j++) {
1034*aad13602SShrirang Abhyankar       /* find row ownership of */
1035*aad13602SShrirang Abhyankar       while (aj[j] >= ranges[proc+1]) proc++;
1036*aad13602SShrirang Abhyankar       nx_all = rranges[proc+1] - rranges[proc];
1037*aad13602SShrirang Abhyankar       col = aj[j] - ranges[proc] + Jranges[proc] + nx_all + nce_all[proc] + nh_all[proc];
1038*aad13602SShrirang Abhyankar       ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr);
1039*aad13602SShrirang Abhyankar     }
1040*aad13602SShrirang Abhyankar     ierr = MatRestoreRow(Jci_xb_trans,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr);
1041*aad13602SShrirang Abhyankar   }
1042*aad13602SShrirang Abhyankar 
1043*aad13602SShrirang Abhyankar   /* 2nd Row block of KKT matrix: [grad Ce, 0, 0, 0] */
1044*aad13602SShrirang Abhyankar   if(pdipm->Ng) {
1045*aad13602SShrirang Abhyankar     ierr = MatGetOwnershipRange(tao->jacobian_equality,&rjstart,NULL);CHKERRQ(ierr);
1046*aad13602SShrirang Abhyankar     for (i=0; i < pdipm->ng; i++){
1047*aad13602SShrirang Abhyankar       row = rstart + pdipm->off_lambdae + i;
1048*aad13602SShrirang Abhyankar 
1049*aad13602SShrirang Abhyankar       ierr = MatGetRow(tao->jacobian_equality,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr);
1050*aad13602SShrirang Abhyankar       proc = 0;
1051*aad13602SShrirang Abhyankar       for (j=0; j < nc; j++) {
1052*aad13602SShrirang Abhyankar         while (aj[j] >= cranges[proc+1]) proc++;
1053*aad13602SShrirang Abhyankar         col = aj[j] - cranges[proc] + Jranges[proc];
1054*aad13602SShrirang Abhyankar         ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr); /* grad g */
1055*aad13602SShrirang Abhyankar       }
1056*aad13602SShrirang Abhyankar       ierr = MatRestoreRow(tao->jacobian_equality,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr);
1057*aad13602SShrirang Abhyankar     }
1058*aad13602SShrirang Abhyankar   }
1059*aad13602SShrirang Abhyankar   /* Jce_xfixed */
1060*aad13602SShrirang Abhyankar   if (pdipm->Nxfixed) {
1061*aad13602SShrirang Abhyankar     ierr = MatGetOwnershipRange(pdipm->Jce_xfixed,&Jcrstart,NULL);CHKERRQ(ierr);
1062*aad13602SShrirang Abhyankar     for (i=0; i < (pdipm->nce - pdipm->ng); i++ ){
1063*aad13602SShrirang Abhyankar       row = rstart + pdipm->off_lambdae + pdipm->ng + i;
1064*aad13602SShrirang Abhyankar 
1065*aad13602SShrirang Abhyankar       ierr = MatGetRow(pdipm->Jce_xfixed,i+Jcrstart,&nc,&cols,NULL);CHKERRQ(ierr);
1066*aad13602SShrirang Abhyankar       if (nc != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"nc != 1");
1067*aad13602SShrirang Abhyankar 
1068*aad13602SShrirang Abhyankar       proc = 0;
1069*aad13602SShrirang Abhyankar       j    = 0;
1070*aad13602SShrirang Abhyankar       while (cols[j] >= cranges[proc+1]) proc++;
1071*aad13602SShrirang Abhyankar       col = cols[j] - cranges[proc] + Jranges[proc];
1072*aad13602SShrirang Abhyankar       ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr);
1073*aad13602SShrirang Abhyankar       ierr = MatRestoreRow(pdipm->Jce_xfixed,i+Jcrstart,&nc,&cols,NULL);CHKERRQ(ierr);
1074*aad13602SShrirang Abhyankar     }
1075*aad13602SShrirang Abhyankar   }
1076*aad13602SShrirang Abhyankar 
1077*aad13602SShrirang Abhyankar   /* 3rd Row block of KKT matrix: [ gradCi, 0, 0, -I] */
1078*aad13602SShrirang Abhyankar   if(pdipm->Nh) {
1079*aad13602SShrirang Abhyankar     ierr = MatGetOwnershipRange(tao->jacobian_inequality,&rjstart,NULL);CHKERRQ(ierr);
1080*aad13602SShrirang Abhyankar     for (i=0; i < pdipm->nh; i++){
1081*aad13602SShrirang Abhyankar       row = rstart + pdipm->off_lambdai + i;
1082*aad13602SShrirang Abhyankar 
1083*aad13602SShrirang Abhyankar       ierr = MatGetRow(tao->jacobian_inequality,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr);
1084*aad13602SShrirang Abhyankar       proc = 0;
1085*aad13602SShrirang Abhyankar       for (j=0; j < nc; j++) {
1086*aad13602SShrirang Abhyankar         while (aj[j] >= cranges[proc+1]) proc++;
1087*aad13602SShrirang Abhyankar         col = aj[j] - cranges[proc] + Jranges[proc];
1088*aad13602SShrirang Abhyankar         ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr); /* grad h */
1089*aad13602SShrirang Abhyankar       }
1090*aad13602SShrirang Abhyankar       ierr = MatRestoreRow(tao->jacobian_inequality,i+rjstart,&nc,&aj,NULL);CHKERRQ(ierr);
1091*aad13602SShrirang Abhyankar     }
1092*aad13602SShrirang Abhyankar     /* -I */
1093*aad13602SShrirang Abhyankar     for (i=0; i < pdipm->nh; i++){
1094*aad13602SShrirang Abhyankar       row = rstart + pdipm->off_lambdai + i;
1095*aad13602SShrirang Abhyankar       col = rstart + pdipm->off_z + i;
1096*aad13602SShrirang Abhyankar       ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr);
1097*aad13602SShrirang Abhyankar     }
1098*aad13602SShrirang Abhyankar   }
1099*aad13602SShrirang Abhyankar 
1100*aad13602SShrirang Abhyankar   /* Jci_xb */
1101*aad13602SShrirang Abhyankar   ierr = MatGetOwnershipRange(pdipm->Jci_xb,&Jcrstart,NULL);CHKERRQ(ierr);
1102*aad13602SShrirang Abhyankar   for (i=0; i < (pdipm->nci - pdipm->nh); i++ ){
1103*aad13602SShrirang Abhyankar     row = rstart + pdipm->off_lambdai + pdipm->nh + i;
1104*aad13602SShrirang Abhyankar 
1105*aad13602SShrirang Abhyankar     ierr = MatGetRow(pdipm->Jci_xb,i+Jcrstart,&nc,&cols,NULL);CHKERRQ(ierr);
1106*aad13602SShrirang Abhyankar     if (nc != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"nc != 1");
1107*aad13602SShrirang Abhyankar     proc = 0;
1108*aad13602SShrirang Abhyankar     for (j=0; j < nc; j++) {
1109*aad13602SShrirang Abhyankar       while (cols[j] >= cranges[proc+1]) proc++;
1110*aad13602SShrirang Abhyankar       col = cols[j] - cranges[proc] + Jranges[proc];
1111*aad13602SShrirang Abhyankar       ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr);
1112*aad13602SShrirang Abhyankar     }
1113*aad13602SShrirang Abhyankar     ierr = MatRestoreRow(pdipm->Jci_xb,i+Jcrstart,&nc,&cols,NULL);CHKERRQ(ierr);
1114*aad13602SShrirang Abhyankar     /* -I */
1115*aad13602SShrirang Abhyankar     col = rstart + pdipm->off_z + pdipm->nh + i;
1116*aad13602SShrirang Abhyankar     ierr = MatPreallocateSet(row,1,&col,dnz,onz);CHKERRQ(ierr);
1117*aad13602SShrirang Abhyankar   }
1118*aad13602SShrirang Abhyankar 
1119*aad13602SShrirang Abhyankar   /* 4-th Row block of KKT matrix: Z and Ci */
1120*aad13602SShrirang Abhyankar   for (i=0; i < pdipm->nci; i++) {
1121*aad13602SShrirang Abhyankar     row     = rstart + pdipm->off_z + i;
1122*aad13602SShrirang Abhyankar     cols1[0] = rstart + pdipm->off_lambdai + i;
1123*aad13602SShrirang Abhyankar     cols1[1] = row;
1124*aad13602SShrirang Abhyankar     ierr = MatPreallocateSet(row,2,cols1,dnz,onz);CHKERRQ(ierr);
1125*aad13602SShrirang Abhyankar   }
1126*aad13602SShrirang Abhyankar 
1127*aad13602SShrirang Abhyankar   /* diagonal entry */
1128*aad13602SShrirang Abhyankar   for (i=0; i<pdipm->n; i++) dnz[i]++; /* diagonal entry */
1129*aad13602SShrirang Abhyankar 
1130*aad13602SShrirang Abhyankar   /* Create KKT matrix */
1131*aad13602SShrirang Abhyankar   ierr = MatCreate(comm,&J);CHKERRQ(ierr);
1132*aad13602SShrirang Abhyankar   ierr = MatSetSizes(J,pdipm->n,pdipm->n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
1133*aad13602SShrirang Abhyankar   ierr = MatSetFromOptions(J);CHKERRQ(ierr);
1134*aad13602SShrirang Abhyankar   ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr);
1135*aad13602SShrirang Abhyankar   ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr);
1136*aad13602SShrirang Abhyankar   /* ierr = MatSetOption(J,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr); */
1137*aad13602SShrirang Abhyankar   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1138*aad13602SShrirang Abhyankar   pdipm->K = J;
1139*aad13602SShrirang Abhyankar 
1140*aad13602SShrirang Abhyankar   /* (8) Set up nonlinear solver SNES */
1141*aad13602SShrirang Abhyankar   ierr = SNESSetFunction(pdipm->snes,NULL,TaoSNESFunction_PDIPM,(void*)tao);CHKERRQ(ierr);
1142*aad13602SShrirang Abhyankar   ierr = SNESSetJacobian(pdipm->snes,J,J,TaoSNESJacobian_PDIPM,(void*)tao);CHKERRQ(ierr);
1143*aad13602SShrirang Abhyankar 
1144*aad13602SShrirang Abhyankar   if (pdipm->solve_reduced_kkt) {
1145*aad13602SShrirang Abhyankar     PC pc;
1146*aad13602SShrirang Abhyankar     ierr = KSPGetPC(tao->ksp,&pc);CHKERRQ(ierr);
1147*aad13602SShrirang Abhyankar     ierr = PCSetType(pc,PCFIELDSPLIT);CHKERRQ(ierr);
1148*aad13602SShrirang Abhyankar     ierr = PCFieldSplitSetType(pc,PC_COMPOSITE_SCHUR);CHKERRQ(ierr);
1149*aad13602SShrirang Abhyankar     ierr = PCFieldSplitSetIS(pc,"2",pdipm->is2);CHKERRQ(ierr);
1150*aad13602SShrirang Abhyankar     ierr = PCFieldSplitSetIS(pc,"1",pdipm->is1);CHKERRQ(ierr);
1151*aad13602SShrirang Abhyankar   }
1152*aad13602SShrirang Abhyankar   ierr = SNESSetFromOptions(pdipm->snes);CHKERRQ(ierr);
1153*aad13602SShrirang Abhyankar 
1154*aad13602SShrirang Abhyankar   /* (9) Insert constant entries to  K */
1155*aad13602SShrirang Abhyankar   /* Set 0.0 to diagonal of K, so that the solver does not complain *about missing diagonal value */
1156*aad13602SShrirang Abhyankar   ierr = MatGetOwnershipRange(J,&rstart,&rend);CHKERRQ(ierr);
1157*aad13602SShrirang Abhyankar   for (i=rstart; i<rend; i++){
1158*aad13602SShrirang Abhyankar     ierr = MatSetValue(J,i,i,0.0,INSERT_VALUES);CHKERRQ(ierr);
1159*aad13602SShrirang Abhyankar   }
1160*aad13602SShrirang Abhyankar 
1161*aad13602SShrirang Abhyankar   /* Row block of K: [ grad Ce, 0, 0, 0] */
1162*aad13602SShrirang Abhyankar   if (pdipm->Nxfixed) {
1163*aad13602SShrirang Abhyankar     ierr = MatGetOwnershipRange(pdipm->Jce_xfixed,&Jcrstart,NULL);CHKERRQ(ierr);
1164*aad13602SShrirang Abhyankar     for (i=0; i < (pdipm->nce - pdipm->ng); i++ ){
1165*aad13602SShrirang Abhyankar       row = rstart + pdipm->off_lambdae + pdipm->ng + i;
1166*aad13602SShrirang Abhyankar 
1167*aad13602SShrirang Abhyankar       ierr = MatGetRow(pdipm->Jce_xfixed,i+Jcrstart,&nc,&cols,&aa);CHKERRQ(ierr);
1168*aad13602SShrirang Abhyankar       proc = 0;
1169*aad13602SShrirang Abhyankar       for (j=0; j < nc; j++) {
1170*aad13602SShrirang Abhyankar         while (cols[j] >= cranges[proc+1]) proc++;
1171*aad13602SShrirang Abhyankar         col = cols[j] - cranges[proc] + Jranges[proc];
1172*aad13602SShrirang Abhyankar         ierr = MatSetValue(J,row,col,aa[j],INSERT_VALUES);CHKERRQ(ierr); /* grad Ce */
1173*aad13602SShrirang Abhyankar         ierr = MatSetValue(J,col,row,aa[j],INSERT_VALUES);CHKERRQ(ierr); /* grad Ce' */
1174*aad13602SShrirang Abhyankar       }
1175*aad13602SShrirang Abhyankar       ierr = MatRestoreRow(pdipm->Jce_xfixed,i+Jcrstart,&nc,&cols,&aa);CHKERRQ(ierr);
1176*aad13602SShrirang Abhyankar     }
1177*aad13602SShrirang Abhyankar   }
1178*aad13602SShrirang Abhyankar 
1179*aad13602SShrirang Abhyankar   /* Row block of K: [ grad Ci, 0, 0, -I] */
1180*aad13602SShrirang Abhyankar   ierr = MatGetOwnershipRange(pdipm->Jci_xb,&Jcrstart,NULL);CHKERRQ(ierr);
1181*aad13602SShrirang Abhyankar   for (i=0; i < (pdipm->nci - pdipm->nh); i++ ){
1182*aad13602SShrirang Abhyankar     row = rstart + pdipm->off_lambdai + pdipm->nh + i;
1183*aad13602SShrirang Abhyankar 
1184*aad13602SShrirang Abhyankar     ierr = MatGetRow(pdipm->Jci_xb,i+Jcrstart,&nc,&cols,&aa);CHKERRQ(ierr);
1185*aad13602SShrirang Abhyankar     proc = 0;
1186*aad13602SShrirang Abhyankar     for (j=0; j < nc; j++) {
1187*aad13602SShrirang Abhyankar       while (cols[j] >= cranges[proc+1]) proc++;
1188*aad13602SShrirang Abhyankar       col = cols[j] - cranges[proc] + Jranges[proc];
1189*aad13602SShrirang Abhyankar       ierr = MatSetValue(J,col,row,-aa[j],INSERT_VALUES);CHKERRQ(ierr);
1190*aad13602SShrirang Abhyankar       ierr = MatSetValue(J,row,col,aa[j],INSERT_VALUES);CHKERRQ(ierr);
1191*aad13602SShrirang Abhyankar     }
1192*aad13602SShrirang Abhyankar     ierr = MatRestoreRow(pdipm->Jci_xb,i+Jcrstart,&nc,&cols,&aa);CHKERRQ(ierr);
1193*aad13602SShrirang Abhyankar 
1194*aad13602SShrirang Abhyankar     col = rstart + pdipm->off_z + pdipm->nh + i;
1195*aad13602SShrirang Abhyankar     ierr = MatSetValue(J,row,col,-1,INSERT_VALUES);CHKERRQ(ierr);
1196*aad13602SShrirang Abhyankar   }
1197*aad13602SShrirang Abhyankar 
1198*aad13602SShrirang Abhyankar   for (i=0; i < pdipm->nh; i++){
1199*aad13602SShrirang Abhyankar     row = rstart + pdipm->off_lambdai + i;
1200*aad13602SShrirang Abhyankar     col = rstart + pdipm->off_z + i;
1201*aad13602SShrirang Abhyankar     ierr = MatSetValue(J,row,col,-1,INSERT_VALUES);CHKERRQ(ierr);
1202*aad13602SShrirang Abhyankar   }
1203*aad13602SShrirang Abhyankar 
1204*aad13602SShrirang Abhyankar   if (pdipm->Nxfixed) {
1205*aad13602SShrirang Abhyankar     ierr = MatDestroy(&Jce_xfixed_trans);CHKERRQ(ierr);
1206*aad13602SShrirang Abhyankar   }
1207*aad13602SShrirang Abhyankar   ierr = MatDestroy(&Jci_xb_trans);CHKERRQ(ierr);
1208*aad13602SShrirang Abhyankar   ierr = PetscFree3(ng_all,nh_all,Jranges);CHKERRQ(ierr);
1209*aad13602SShrirang Abhyankar   PetscFunctionReturn(0);
1210*aad13602SShrirang Abhyankar }
1211*aad13602SShrirang Abhyankar 
1212*aad13602SShrirang Abhyankar /*
1213*aad13602SShrirang Abhyankar    TaoDestroy_PDIPM - Destroys the pdipm object
1214*aad13602SShrirang Abhyankar 
1215*aad13602SShrirang Abhyankar    Input:
1216*aad13602SShrirang Abhyankar    full pdipm
1217*aad13602SShrirang Abhyankar 
1218*aad13602SShrirang Abhyankar    Output:
1219*aad13602SShrirang Abhyankar    Destroyed pdipm
1220*aad13602SShrirang Abhyankar */
1221*aad13602SShrirang Abhyankar PetscErrorCode TaoDestroy_PDIPM(Tao tao)
1222*aad13602SShrirang Abhyankar {
1223*aad13602SShrirang Abhyankar   TAO_PDIPM      *pdipm = (TAO_PDIPM*)tao->data;
1224*aad13602SShrirang Abhyankar   PetscErrorCode ierr;
1225*aad13602SShrirang Abhyankar 
1226*aad13602SShrirang Abhyankar   PetscFunctionBegin;
1227*aad13602SShrirang Abhyankar   /* Freeing Vectors assocaiated with KKT (X) */
1228*aad13602SShrirang Abhyankar   ierr = VecDestroy(&pdipm->x);CHKERRQ(ierr); /* Solution x */
1229*aad13602SShrirang Abhyankar   ierr = VecDestroy(&pdipm->lambdae);CHKERRQ(ierr); /* Equality constraints lagrangian multiplier*/
1230*aad13602SShrirang Abhyankar   ierr = VecDestroy(&pdipm->lambdai);CHKERRQ(ierr); /* Inequality constraints lagrangian multiplier*/
1231*aad13602SShrirang Abhyankar   ierr = VecDestroy(&pdipm->z);CHKERRQ(ierr);       /* Slack variables */
1232*aad13602SShrirang Abhyankar   ierr = VecDestroy(&pdipm->X);CHKERRQ(ierr);       /* Big KKT system vector [x; lambdae; lambdai; z] */
1233*aad13602SShrirang Abhyankar 
1234*aad13602SShrirang Abhyankar   /* work vectors */
1235*aad13602SShrirang Abhyankar   ierr = VecDestroy(&pdipm->lambdae_xfixed);CHKERRQ(ierr);
1236*aad13602SShrirang Abhyankar   ierr = VecDestroy(&pdipm->lambdai_xb);CHKERRQ(ierr);
1237*aad13602SShrirang Abhyankar 
1238*aad13602SShrirang Abhyankar   /* Legrangian equality and inequality Vec */
1239*aad13602SShrirang Abhyankar   ierr = VecDestroy(&pdipm->ce);CHKERRQ(ierr); /* Vec of equality constraints */
1240*aad13602SShrirang Abhyankar   ierr = VecDestroy(&pdipm->ci);CHKERRQ(ierr); /* Vec of inequality constraints */
1241*aad13602SShrirang Abhyankar 
1242*aad13602SShrirang Abhyankar   /* Matrices */
1243*aad13602SShrirang Abhyankar   ierr = MatDestroy(&pdipm->Jce_xfixed);CHKERRQ(ierr);
1244*aad13602SShrirang Abhyankar   ierr = MatDestroy(&pdipm->Jci_xb);CHKERRQ(ierr); /* Jacobian of inequality constraints Jci = [tao->jacobian_inequality ; J(nxub); J(nxlb); J(nxbx)] */
1245*aad13602SShrirang Abhyankar   ierr = MatDestroy(&pdipm->K);CHKERRQ(ierr);
1246*aad13602SShrirang Abhyankar 
1247*aad13602SShrirang Abhyankar   /* Index Sets */
1248*aad13602SShrirang Abhyankar   if (pdipm->Nxub) {
1249*aad13602SShrirang Abhyankar     ierr = ISDestroy(&pdipm->isxub);CHKERRQ(ierr);    /* Finite upper bound only -inf < x < ub */
1250*aad13602SShrirang Abhyankar   }
1251*aad13602SShrirang Abhyankar 
1252*aad13602SShrirang Abhyankar   if (pdipm->Nxlb) {
1253*aad13602SShrirang Abhyankar     ierr = ISDestroy(&pdipm->isxlb);CHKERRQ(ierr);    /* Finite lower bound only  lb <= x < inf */
1254*aad13602SShrirang Abhyankar   }
1255*aad13602SShrirang Abhyankar 
1256*aad13602SShrirang Abhyankar   if (pdipm->Nxfixed) {
1257*aad13602SShrirang Abhyankar     ierr = ISDestroy(&pdipm->isxfixed);CHKERRQ(ierr); /* Fixed variables         lb =  x = ub */
1258*aad13602SShrirang Abhyankar   }
1259*aad13602SShrirang Abhyankar 
1260*aad13602SShrirang Abhyankar   if (pdipm->Nxbox) {
1261*aad13602SShrirang Abhyankar     ierr = ISDestroy(&pdipm->isxbox);CHKERRQ(ierr);   /* Boxed variables         lb <= x <= ub */
1262*aad13602SShrirang Abhyankar   }
1263*aad13602SShrirang Abhyankar 
1264*aad13602SShrirang Abhyankar   if (pdipm->Nxfree) {
1265*aad13602SShrirang Abhyankar     ierr = ISDestroy(&pdipm->isxfree);CHKERRQ(ierr);  /* Free variables        -inf <= x <= inf */
1266*aad13602SShrirang Abhyankar   }
1267*aad13602SShrirang Abhyankar 
1268*aad13602SShrirang Abhyankar   if (pdipm->solve_reduced_kkt) {
1269*aad13602SShrirang Abhyankar     ierr = ISDestroy(&pdipm->is1);CHKERRQ(ierr);
1270*aad13602SShrirang Abhyankar     ierr = ISDestroy(&pdipm->is2);CHKERRQ(ierr);
1271*aad13602SShrirang Abhyankar   }
1272*aad13602SShrirang Abhyankar 
1273*aad13602SShrirang Abhyankar   /* SNES */
1274*aad13602SShrirang Abhyankar   ierr = SNESDestroy(&pdipm->snes);CHKERRQ(ierr); /* Nonlinear solver */
1275*aad13602SShrirang Abhyankar   ierr = PetscFree(pdipm->nce_all);CHKERRQ(ierr);
1276*aad13602SShrirang Abhyankar   ierr = MatDestroy(&pdipm->jac_equality_trans);CHKERRQ(ierr);
1277*aad13602SShrirang Abhyankar   ierr = MatDestroy(&pdipm->jac_inequality_trans);CHKERRQ(ierr);
1278*aad13602SShrirang Abhyankar 
1279*aad13602SShrirang Abhyankar   /* Destroy pdipm */
1280*aad13602SShrirang Abhyankar   ierr = PetscFree(tao->data);CHKERRQ(ierr); /* Holding locations of pdipm */
1281*aad13602SShrirang Abhyankar 
1282*aad13602SShrirang Abhyankar   /* Destroy Dual */
1283*aad13602SShrirang Abhyankar   ierr = VecDestroy(&tao->DE);CHKERRQ(ierr); /* equality dual */
1284*aad13602SShrirang Abhyankar   ierr = VecDestroy(&tao->DI);CHKERRQ(ierr); /* dinequality dual */
1285*aad13602SShrirang Abhyankar   PetscFunctionReturn(0);
1286*aad13602SShrirang Abhyankar }
1287*aad13602SShrirang Abhyankar 
1288*aad13602SShrirang Abhyankar PetscErrorCode TaoSetFromOptions_PDIPM(PetscOptionItems *PetscOptionsObject,Tao tao)
1289*aad13602SShrirang Abhyankar {
1290*aad13602SShrirang Abhyankar   TAO_PDIPM      *pdipm = (TAO_PDIPM*)tao->data;
1291*aad13602SShrirang Abhyankar   PetscErrorCode ierr;
1292*aad13602SShrirang Abhyankar 
1293*aad13602SShrirang Abhyankar   PetscFunctionBegin;
1294*aad13602SShrirang Abhyankar   ierr = PetscOptionsHead(PetscOptionsObject,"PDIPM method for constrained optimization");CHKERRQ(ierr);
1295*aad13602SShrirang Abhyankar   ierr = PetscOptionsReal("-tao_pdipm_push_init_slack","parameter to push initial slack variables away from bounds",NULL,pdipm->push_init_slack,&pdipm->push_init_slack,NULL);CHKERRQ(ierr);
1296*aad13602SShrirang Abhyankar   ierr = PetscOptionsReal("-tao_pdipm_push_init_lambdai","parameter to push initial (inequality) dual variables away from bounds",NULL,pdipm->push_init_lambdai,&pdipm->push_init_lambdai,NULL);CHKERRQ(ierr);
1297*aad13602SShrirang Abhyankar   ierr = PetscOptionsBool("-tao_pdipm_solve_reduced_kkt","Solve reduced KKT system using Schur-complement",NULL,pdipm->solve_reduced_kkt,&pdipm->solve_reduced_kkt,NULL);CHKERRQ(ierr);
1298*aad13602SShrirang Abhyankar   ierr = PetscOptionsReal("-tao_pdipm_mu_update_factor","Update scalar for barrier parameter (mu) update",NULL,pdipm->mu_update_factor,&pdipm->mu_update_factor,NULL);CHKERRQ(ierr);
1299*aad13602SShrirang Abhyankar   ierr = PetscOptionsTail();CHKERRQ(ierr);
1300*aad13602SShrirang Abhyankar   PetscFunctionReturn(0);
1301*aad13602SShrirang Abhyankar }
1302*aad13602SShrirang Abhyankar 
1303*aad13602SShrirang Abhyankar /*MC
1304*aad13602SShrirang Abhyankar   TAOPDIPM - Barrier-based primal-dual interior point algorithm for generally constrained optimization.
1305*aad13602SShrirang Abhyankar 
1306*aad13602SShrirang Abhyankar   Option Database Keys:
1307*aad13602SShrirang Abhyankar +   -tao_pdipm_push_init_lambdai - parameter to push initial dual variables away from bounds (> 0)
1308*aad13602SShrirang Abhyankar .   -tao_pdipm_push_init_slack  - parameter to push initial slack variables away from bounds (> 0)
1309*aad13602SShrirang Abhyankar -   -tao_pdipm_mu_update_factor - update scalar for barrier parameter (mu) update (> 0)
1310*aad13602SShrirang Abhyankar 
1311*aad13602SShrirang Abhyankar   Level: beginner
1312*aad13602SShrirang Abhyankar M*/
1313*aad13602SShrirang Abhyankar PETSC_EXTERN PetscErrorCode TaoCreate_PDIPM(Tao tao)
1314*aad13602SShrirang Abhyankar {
1315*aad13602SShrirang Abhyankar   TAO_PDIPM      *pdipm;
1316*aad13602SShrirang Abhyankar   PetscErrorCode ierr;
1317*aad13602SShrirang Abhyankar 
1318*aad13602SShrirang Abhyankar   PetscFunctionBegin;
1319*aad13602SShrirang Abhyankar   tao->ops->setup          = TaoSetup_PDIPM;
1320*aad13602SShrirang Abhyankar   tao->ops->solve          = TaoSolve_PDIPM;
1321*aad13602SShrirang Abhyankar   tao->ops->setfromoptions = TaoSetFromOptions_PDIPM;
1322*aad13602SShrirang Abhyankar   tao->ops->destroy        = TaoDestroy_PDIPM;
1323*aad13602SShrirang Abhyankar 
1324*aad13602SShrirang Abhyankar   ierr = PetscNewLog(tao,&pdipm);CHKERRQ(ierr);
1325*aad13602SShrirang Abhyankar   tao->data = (void*)pdipm;
1326*aad13602SShrirang Abhyankar 
1327*aad13602SShrirang Abhyankar   pdipm->nx      = pdipm->Nx      = 0;
1328*aad13602SShrirang Abhyankar   pdipm->nxfixed = pdipm->Nxfixed = 0;
1329*aad13602SShrirang Abhyankar   pdipm->nxlb    = pdipm->Nxlb    = 0;
1330*aad13602SShrirang Abhyankar   pdipm->nxub    = pdipm->Nxub    = 0;
1331*aad13602SShrirang Abhyankar   pdipm->nxbox   = pdipm->Nxbox   = 0;
1332*aad13602SShrirang Abhyankar   pdipm->nxfree  = pdipm->Nxfree  = 0;
1333*aad13602SShrirang Abhyankar 
1334*aad13602SShrirang Abhyankar   pdipm->ng = pdipm->Ng = pdipm->nce = pdipm->Nce = 0;
1335*aad13602SShrirang Abhyankar   pdipm->nh = pdipm->Nh = pdipm->nci = pdipm->Nci = 0;
1336*aad13602SShrirang Abhyankar   pdipm->n  = pdipm->N  = 0;
1337*aad13602SShrirang Abhyankar   pdipm->mu = 1.0;
1338*aad13602SShrirang Abhyankar   pdipm->mu_update_factor = 0.1;
1339*aad13602SShrirang Abhyankar 
1340*aad13602SShrirang Abhyankar   pdipm->push_init_slack   = 1.0;
1341*aad13602SShrirang Abhyankar   pdipm->push_init_lambdai = 1.0;
1342*aad13602SShrirang Abhyankar   pdipm->solve_reduced_kkt = PETSC_FALSE;
1343*aad13602SShrirang Abhyankar 
1344*aad13602SShrirang Abhyankar   /* Override default settings (unless already changed) */
1345*aad13602SShrirang Abhyankar   if (!tao->max_it_changed) tao->max_it = 200;
1346*aad13602SShrirang Abhyankar   if (!tao->max_funcs_changed) tao->max_funcs = 500;
1347*aad13602SShrirang Abhyankar 
1348*aad13602SShrirang Abhyankar   ierr = SNESCreate(((PetscObject)tao)->comm,&pdipm->snes);CHKERRQ(ierr);
1349*aad13602SShrirang Abhyankar   ierr = SNESSetOptionsPrefix(pdipm->snes,tao->hdr.prefix);CHKERRQ(ierr);
1350*aad13602SShrirang Abhyankar   ierr = SNESGetKSP(pdipm->snes,&tao->ksp);CHKERRQ(ierr);
1351*aad13602SShrirang Abhyankar   ierr = PetscObjectReference((PetscObject)tao->ksp);CHKERRQ(ierr);
1352*aad13602SShrirang Abhyankar   PetscFunctionReturn(0);
1353*aad13602SShrirang Abhyankar }
1354