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