xref: /petsc/src/tao/constrained/tutorials/ex1.c (revision 661095bbfddda9a1493a32ea0d2305a96eb189ff)
1aad13602SShrirang Abhyankar /* Program usage: mpiexec -n 2 ./ex1 [-help] [all TAO options] */
2aad13602SShrirang Abhyankar 
3aad13602SShrirang Abhyankar /* ----------------------------------------------------------------------
4aad13602SShrirang Abhyankar min f = (x0 - 2)^2 + (x1 - 2)^2 - 2*(x0 + x1)
5aad13602SShrirang Abhyankar s.t.  x0^2 + x1 - 2 = 0
6aad13602SShrirang Abhyankar       0  <= x0^2 - x1 <= 1
7aad13602SShrirang Abhyankar       -1 <= x0, x1 <= 2
8aad13602SShrirang Abhyankar ---------------------------------------------------------------------- */
9aad13602SShrirang Abhyankar 
10aad13602SShrirang Abhyankar #include <petsctao.h>
11aad13602SShrirang Abhyankar 
12aad13602SShrirang Abhyankar static  char help[]= "Solves constrained optimiztion problem using pdipm.\n\
13aad13602SShrirang Abhyankar Input parameters include:\n\
14aad13602SShrirang Abhyankar   -tao_type pdipm    : sets Tao solver\n\
158e58fa1dSresundermann   -no_eq             : removes the equaility constraints from the problem\n\
16aad13602SShrirang Abhyankar   -snes_fd           : snes with finite difference Jacobian (needed for pdipm)\n\
17aad13602SShrirang Abhyankar   -tao_cmonitor      : convergence monitor with constraint norm \n\
18aad13602SShrirang Abhyankar   -tao_view_solution : view exact solution at each itteration\n\
1912d688e0SRylee Sundermann   Note: external package mumps is requried to run either for pdipm. Additionally This is designed for a maximum of 2 processors, the code will error if size > 2.\n\n";
20aad13602SShrirang Abhyankar 
21aad13602SShrirang Abhyankar /*
22aad13602SShrirang Abhyankar    User-defined application context - contains data needed by the
23aad13602SShrirang Abhyankar    application-provided call-back routines, FormFunction(),
24aad13602SShrirang Abhyankar    FormGradient(), and FormHessian().
25aad13602SShrirang Abhyankar */
26aad13602SShrirang Abhyankar 
27aad13602SShrirang Abhyankar /*
28aad13602SShrirang Abhyankar    x,d in R^n
29aad13602SShrirang Abhyankar    f in R
30aad13602SShrirang Abhyankar    bin in R^mi
31aad13602SShrirang Abhyankar    beq in R^me
32aad13602SShrirang Abhyankar    Aeq in R^(me x n)
33aad13602SShrirang Abhyankar    Ain in R^(mi x n)
34aad13602SShrirang Abhyankar    H in R^(n x n)
35aad13602SShrirang Abhyankar    min f=(1/2)*x'*H*x + d'*x
36aad13602SShrirang Abhyankar    s.t.  Aeq*x == beq
37aad13602SShrirang Abhyankar          Ain*x >= bin
38aad13602SShrirang Abhyankar */
39aad13602SShrirang Abhyankar typedef struct {
40aad13602SShrirang Abhyankar   PetscInt   n;  /* Global length of x */
41aad13602SShrirang Abhyankar   PetscInt   ne; /* Global number of equality constraints */
42aad13602SShrirang Abhyankar   PetscInt   ni; /* Global number of inequality constraints */
438e58fa1dSresundermann   PetscBool  noeqflag;
44aad13602SShrirang Abhyankar   Vec        x,xl,xu;
45aad13602SShrirang Abhyankar   Vec        ce,ci,bl,bu,Xseq;
46aad13602SShrirang Abhyankar   Mat        Ae,Ai,H;
47aad13602SShrirang Abhyankar   VecScatter scat;
48aad13602SShrirang Abhyankar } AppCtx;
49aad13602SShrirang Abhyankar 
50*661095bbSAlp Dener 
51aad13602SShrirang Abhyankar /* -------- User-defined Routines --------- */
52aad13602SShrirang Abhyankar PetscErrorCode InitializeProblem(AppCtx *);
53aad13602SShrirang Abhyankar PetscErrorCode DestroyProblem(AppCtx *);
54aad13602SShrirang Abhyankar PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal *,Vec,void *);
55aad13602SShrirang Abhyankar PetscErrorCode FormHessian(Tao,Vec,Mat,Mat, void*);
56aad13602SShrirang Abhyankar PetscErrorCode FormInequalityConstraints(Tao,Vec,Vec,void*);
57aad13602SShrirang Abhyankar PetscErrorCode FormEqualityConstraints(Tao,Vec,Vec,void*);
58aad13602SShrirang Abhyankar PetscErrorCode FormInequalityJacobian(Tao,Vec,Mat,Mat, void*);
59aad13602SShrirang Abhyankar PetscErrorCode FormEqualityJacobian(Tao,Vec,Mat,Mat, void*);
60aad13602SShrirang Abhyankar 
61aad13602SShrirang Abhyankar PetscErrorCode main(int argc,char **argv)
62aad13602SShrirang Abhyankar {
63aad13602SShrirang Abhyankar   PetscErrorCode ierr;  /* used to check for functions returning nonzeros */
64aad13602SShrirang Abhyankar   Tao            tao;
65aad13602SShrirang Abhyankar   KSP            ksp;
66aad13602SShrirang Abhyankar   PC             pc;
67aad13602SShrirang Abhyankar   AppCtx         user;  /* application context */
68aad13602SShrirang Abhyankar   Vec            x;
69aad13602SShrirang Abhyankar   PetscMPIInt    rank;
70*661095bbSAlp Dener   TaoType        type;
71*661095bbSAlp Dener   PetscBool      pdipm;
72aad13602SShrirang Abhyankar 
73aad13602SShrirang Abhyankar   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
74ffc4695bSBarry Smith   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr);
75aad13602SShrirang Abhyankar   if (rank>1){
76aad13602SShrirang Abhyankar     SETERRQ(PETSC_COMM_SELF,1,"More than 2 processors detected. Example written to use max of 2 processors.\n");
77aad13602SShrirang Abhyankar   }
78aad13602SShrirang Abhyankar 
79aad13602SShrirang Abhyankar   ierr = PetscPrintf(PETSC_COMM_WORLD,"---- Constrained Problem -----\n");CHKERRQ(ierr);
80aad13602SShrirang Abhyankar   ierr = InitializeProblem(&user);CHKERRQ(ierr); /* sets up problem, function below */
81aad13602SShrirang Abhyankar   ierr = TaoCreate(PETSC_COMM_WORLD,&tao);CHKERRQ(ierr);
82aad13602SShrirang Abhyankar   ierr = TaoSetType(tao,TAOPDIPM);CHKERRQ(ierr);
83aad13602SShrirang Abhyankar   ierr = TaoSetInitialVector(tao,user.x);CHKERRQ(ierr); /* gets solution vector from problem */
84aad13602SShrirang Abhyankar   ierr = TaoSetVariableBounds(tao,user.xl,user.xu);CHKERRQ(ierr); /* sets lower upper bounds from given solution */
85aad13602SShrirang Abhyankar   ierr = TaoSetObjectiveAndGradientRoutine(tao,FormFunctionGradient,(void*)&user);CHKERRQ(ierr);
86aad13602SShrirang Abhyankar 
878e58fa1dSresundermann   if (!user.noeqflag){
88aad13602SShrirang Abhyankar     ierr = TaoSetEqualityConstraintsRoutine(tao,user.ce,FormEqualityConstraints,(void*)&user);CHKERRQ(ierr);
898e58fa1dSresundermann   }
90aad13602SShrirang Abhyankar   ierr = TaoSetInequalityConstraintsRoutine(tao,user.ci,FormInequalityConstraints,(void*)&user);CHKERRQ(ierr);
91aad13602SShrirang Abhyankar 
928e58fa1dSresundermann   if (!user.noeqflag){
93aad13602SShrirang Abhyankar     ierr = TaoSetJacobianEqualityRoutine(tao,user.Ae,user.Ae,FormEqualityJacobian,(void*)&user);CHKERRQ(ierr); /* equality jacobian */
948e58fa1dSresundermann   }
95aad13602SShrirang Abhyankar   ierr = TaoSetJacobianInequalityRoutine(tao,user.Ai,user.Ai,FormInequalityJacobian,(void*)&user);CHKERRQ(ierr); /* inequality jacobian */
96aad13602SShrirang Abhyankar   ierr = TaoSetTolerances(tao,1.e-6,1.e-6,1.e-6);CHKERRQ(ierr);
97aad13602SShrirang Abhyankar   ierr = TaoSetConstraintTolerances(tao,1.e-6,1.e-6);CHKERRQ(ierr);
98aad13602SShrirang Abhyankar 
99aad13602SShrirang Abhyankar   ierr = TaoGetKSP(tao,&ksp);CHKERRQ(ierr);
100aad13602SShrirang Abhyankar   ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
10112d688e0SRylee Sundermann   ierr = PCSetType(pc,PCCHOLESKY);CHKERRQ(ierr);
102aad13602SShrirang Abhyankar   /*
103aad13602SShrirang Abhyankar       This algorithm produces matrices with zeros along the diagonal therefore we use
10412d688e0SRylee Sundermann     MUMPS which provides solver for indefinite matrices
105aad13602SShrirang Abhyankar   */
10612d688e0SRylee Sundermann   ierr = PCFactorSetMatSolverType(pc,MATSOLVERMUMPS);CHKERRQ(ierr);  /* requires mumps to solve pdipm */
107aad13602SShrirang Abhyankar   ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr);
108aad13602SShrirang Abhyankar   ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr);
109aad13602SShrirang Abhyankar 
110aad13602SShrirang Abhyankar   ierr = TaoSetFromOptions(tao);CHKERRQ(ierr);
111*661095bbSAlp Dener   ierr = TaoGetType(tao,&type);CHKERRQ(ierr);
112*661095bbSAlp Dener   ierr = PetscObjectTypeCompare((PetscObject)tao,TAOPDIPM,&pdipm);CHKERRQ(ierr);
113*661095bbSAlp Dener   if (pdipm) {
114*661095bbSAlp Dener     ierr = TaoSetHessianRoutine(tao,user.H,user.H,FormHessian,(void*)&user);CHKERRQ(ierr);
115*661095bbSAlp Dener   }
116aad13602SShrirang Abhyankar 
117aad13602SShrirang Abhyankar   ierr = TaoSolve(tao);CHKERRQ(ierr);
118aad13602SShrirang Abhyankar   ierr = TaoGetSolutionVector(tao,&x);CHKERRQ(ierr);
119aad13602SShrirang Abhyankar   ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
120aad13602SShrirang Abhyankar 
121aad13602SShrirang Abhyankar   /* Free objects */
122aad13602SShrirang Abhyankar   ierr = DestroyProblem(&user);CHKERRQ(ierr);
123aad13602SShrirang Abhyankar   ierr = TaoDestroy(&tao);CHKERRQ(ierr);
124aad13602SShrirang Abhyankar   ierr = PetscFinalize();
125aad13602SShrirang Abhyankar   return ierr;
126aad13602SShrirang Abhyankar }
127aad13602SShrirang Abhyankar 
128aad13602SShrirang Abhyankar PetscErrorCode InitializeProblem(AppCtx *user)
129aad13602SShrirang Abhyankar {
130aad13602SShrirang Abhyankar   PetscErrorCode ierr;
131aad13602SShrirang Abhyankar   PetscMPIInt    size;
132aad13602SShrirang Abhyankar   PetscMPIInt    rank;
133aad13602SShrirang Abhyankar   PetscInt       nloc,neloc,niloc;
134aad13602SShrirang Abhyankar 
135aad13602SShrirang Abhyankar   PetscFunctionBegin;
136ffc4695bSBarry Smith   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
137ffc4695bSBarry Smith   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr);
1388e58fa1dSresundermann   user->noeqflag = PETSC_FALSE;
1398e58fa1dSresundermann   ierr = PetscOptionsGetBool(NULL,NULL,"-no_eq",&user->noeqflag,NULL);CHKERRQ(ierr);
1408e58fa1dSresundermann   if (!user->noeqflag) {
1418e58fa1dSresundermann     ierr = PetscPrintf(PETSC_COMM_WORLD,"Solution should be f(1,1)=-2\n");CHKERRQ(ierr);
1428e58fa1dSresundermann   }
143aad13602SShrirang Abhyankar 
1442d4ee042Sprj-   /* create vector x and set initial values */
145aad13602SShrirang Abhyankar   user->n = 2; /* global length */
146aad13602SShrirang Abhyankar   nloc = (rank==0)?user->n:0;
147*661095bbSAlp Dener   ierr = VecCreate(PETSC_COMM_WORLD,&user->x);CHKERRQ(ierr);
148*661095bbSAlp Dener   ierr = VecSetSizes(user->x,nloc,user->n);CHKERRQ(ierr);
149*661095bbSAlp Dener   ierr = VecSetFromOptions(user->x);CHKERRQ(ierr);
150aad13602SShrirang Abhyankar   ierr = VecSet(user->x,0);CHKERRQ(ierr);
151aad13602SShrirang Abhyankar 
152aad13602SShrirang Abhyankar   /* create and set lower and upper bound vectors */
153aad13602SShrirang Abhyankar   ierr = VecDuplicate(user->x,&user->xl);CHKERRQ(ierr);
154aad13602SShrirang Abhyankar   ierr = VecDuplicate(user->x,&user->xu);CHKERRQ(ierr);
155aad13602SShrirang Abhyankar   ierr = VecSet(user->xl,-1.0);CHKERRQ(ierr);
156aad13602SShrirang Abhyankar   ierr = VecSet(user->xu,2.0);CHKERRQ(ierr);
157aad13602SShrirang Abhyankar 
158aad13602SShrirang Abhyankar   /* create scater to zero */
159aad13602SShrirang Abhyankar   ierr = VecScatterCreateToZero(user->x,&user->scat,&user->Xseq);CHKERRQ(ierr);
160aad13602SShrirang Abhyankar 
161aad13602SShrirang Abhyankar     user->ne = 1;
162aad13602SShrirang Abhyankar     neloc = (rank==0)?user->ne:0;
163aad13602SShrirang Abhyankar 
1648e58fa1dSresundermann   if (!user->noeqflag){
165*661095bbSAlp Dener     ierr = VecCreate(PETSC_COMM_WORLD,&user->ce);CHKERRQ(ierr); /* a 1x1 vec for equality constraints */
166*661095bbSAlp Dener     ierr = VecSetSizes(user->ce,neloc,user->ne);CHKERRQ(ierr);
167*661095bbSAlp Dener     ierr = VecSetFromOptions(user->ce);CHKERRQ(ierr);
168*661095bbSAlp Dener     ierr = VecSetUp(user->ce);CHKERRQ(ierr);
1698e58fa1dSresundermann   }
170aad13602SShrirang Abhyankar   user->ni = 2;
171aad13602SShrirang Abhyankar   niloc = (rank==0)?user->ni:0;
172*661095bbSAlp Dener   ierr = VecCreate(PETSC_COMM_WORLD,&user->ci);CHKERRQ(ierr); /* a 2x1 vec for inequality constraints */
173*661095bbSAlp Dener   ierr = VecSetSizes(user->ci,niloc,user->ni);CHKERRQ(ierr);
174*661095bbSAlp Dener   ierr = VecSetFromOptions(user->ci);CHKERRQ(ierr);
175*661095bbSAlp Dener   ierr = VecSetUp(user->ci);CHKERRQ(ierr);
176aad13602SShrirang Abhyankar 
177aad13602SShrirang Abhyankar   /* nexn & nixn matricies for equaly and inequalty constriants */
1788e58fa1dSresundermann   if (!user->noeqflag){
179aad13602SShrirang Abhyankar     ierr = MatCreate(PETSC_COMM_WORLD,&user->Ae);CHKERRQ(ierr);
1808e58fa1dSresundermann     ierr = MatSetSizes(user->Ae,neloc,nloc,user->ne,user->n);CHKERRQ(ierr);
1818e58fa1dSresundermann     ierr = MatSetFromOptions(user->Ae);CHKERRQ(ierr);
182*661095bbSAlp Dener     ierr = MatSetUp(user->Ae);CHKERRQ(ierr);
1838e58fa1dSresundermann   }
1848e58fa1dSresundermann 
185aad13602SShrirang Abhyankar   ierr = MatCreate(PETSC_COMM_WORLD,&user->Ai);CHKERRQ(ierr);
186aad13602SShrirang Abhyankar   ierr = MatCreate(PETSC_COMM_WORLD,&user->H);CHKERRQ(ierr);
187aad13602SShrirang Abhyankar 
188aad13602SShrirang Abhyankar   ierr = MatSetSizes(user->Ai,niloc,nloc,user->ni,user->n);CHKERRQ(ierr);
189aad13602SShrirang Abhyankar   ierr = MatSetSizes(user->H,nloc,nloc,user->n,user->n);CHKERRQ(ierr);
190aad13602SShrirang Abhyankar 
191aad13602SShrirang Abhyankar   ierr = MatSetFromOptions(user->Ai);CHKERRQ(ierr);
192aad13602SShrirang Abhyankar   ierr = MatSetFromOptions(user->H);CHKERRQ(ierr);
193*661095bbSAlp Dener 
194*661095bbSAlp Dener   ierr = MatSetUp(user->Ai);CHKERRQ(ierr);
195*661095bbSAlp Dener   ierr = MatSetUp(user->H);CHKERRQ(ierr);
196aad13602SShrirang Abhyankar   PetscFunctionReturn(0);
197aad13602SShrirang Abhyankar }
198aad13602SShrirang Abhyankar 
199aad13602SShrirang Abhyankar PetscErrorCode DestroyProblem(AppCtx *user)
200aad13602SShrirang Abhyankar {
201aad13602SShrirang Abhyankar   PetscErrorCode ierr;
202aad13602SShrirang Abhyankar 
203aad13602SShrirang Abhyankar   PetscFunctionBegin;
2048e58fa1dSresundermann   if (!user->noeqflag){
205aad13602SShrirang Abhyankar    ierr = MatDestroy(&user->Ae);CHKERRQ(ierr);
2068e58fa1dSresundermann   }
207aad13602SShrirang Abhyankar   ierr = MatDestroy(&user->Ai);CHKERRQ(ierr);
208aad13602SShrirang Abhyankar   ierr = MatDestroy(&user->H);CHKERRQ(ierr);
209aad13602SShrirang Abhyankar 
210aad13602SShrirang Abhyankar   ierr = VecDestroy(&user->x);CHKERRQ(ierr);
2118e58fa1dSresundermann   if (!user->noeqflag){
212aad13602SShrirang Abhyankar     ierr = VecDestroy(&user->ce);CHKERRQ(ierr);
2138e58fa1dSresundermann   }
214aad13602SShrirang Abhyankar   ierr = VecDestroy(&user->ci);CHKERRQ(ierr);
215aad13602SShrirang Abhyankar   ierr = VecDestroy(&user->xl);CHKERRQ(ierr);
216aad13602SShrirang Abhyankar   ierr = VecDestroy(&user->xu);CHKERRQ(ierr);
217aad13602SShrirang Abhyankar   ierr = VecDestroy(&user->Xseq);CHKERRQ(ierr);
218aad13602SShrirang Abhyankar   ierr = VecScatterDestroy(&user->scat);CHKERRQ(ierr);
219aad13602SShrirang Abhyankar   PetscFunctionReturn(0);
220aad13602SShrirang Abhyankar }
221aad13602SShrirang Abhyankar 
222aad13602SShrirang Abhyankar /*
223aad13602SShrirang Abhyankar   f(X) = (x0 - 2)^2 + (x1 - 2)^2 - 2*(x0 + x1)
224aad13602SShrirang Abhyankar   G(X) = fx = [2.0*(x[0]-2.0) - 2.0; 2.0*(x[1]-2.0) - 2.0]
225aad13602SShrirang Abhyankar */
226aad13602SShrirang Abhyankar PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *f, Vec G, void *ctx)
227aad13602SShrirang Abhyankar {
228aad13602SShrirang Abhyankar   PetscScalar       g;
229aad13602SShrirang Abhyankar   const PetscScalar *x;
230aad13602SShrirang Abhyankar   MPI_Comm          comm;
231aad13602SShrirang Abhyankar   PetscMPIInt       rank;
232aad13602SShrirang Abhyankar   PetscErrorCode    ierr;
233aad13602SShrirang Abhyankar   PetscReal         fin;
234aad13602SShrirang Abhyankar   AppCtx            *user=(AppCtx*)ctx;
235aad13602SShrirang Abhyankar   Vec               Xseq=user->Xseq;
236aad13602SShrirang Abhyankar   VecScatter        scat=user->scat;
237aad13602SShrirang Abhyankar 
238aad13602SShrirang Abhyankar   PetscFunctionBegin;
239aad13602SShrirang Abhyankar   ierr = PetscObjectGetComm((PetscObject)tao,&comm);CHKERRQ(ierr);
240ffc4695bSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr);
241aad13602SShrirang Abhyankar 
242aad13602SShrirang Abhyankar   ierr = VecScatterBegin(scat,X,Xseq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
243aad13602SShrirang Abhyankar   ierr = VecScatterEnd(scat,X,Xseq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
244aad13602SShrirang Abhyankar 
245aad13602SShrirang Abhyankar   fin = 0.0;
246aad13602SShrirang Abhyankar   if (!rank) {
247aad13602SShrirang Abhyankar     ierr = VecGetArrayRead(Xseq,&x);CHKERRQ(ierr);
248aad13602SShrirang Abhyankar     fin = (x[0]-2.0)*(x[0]-2.0) + (x[1]-2.0)*(x[1]-2.0) - 2.0*(x[0]+x[1]);
249aad13602SShrirang Abhyankar     g = 2.0*(x[0]-2.0) - 2.0;
250aad13602SShrirang Abhyankar     ierr = VecSetValue(G,0,g,INSERT_VALUES);CHKERRQ(ierr);
251aad13602SShrirang Abhyankar     g = 2.0*(x[1]-2.0) - 2.0;
252aad13602SShrirang Abhyankar     ierr = VecSetValue(G,1,g,INSERT_VALUES);CHKERRQ(ierr);
253aad13602SShrirang Abhyankar     ierr = VecRestoreArrayRead(Xseq,&x);CHKERRQ(ierr);
254aad13602SShrirang Abhyankar   }
255ffc4695bSBarry Smith   ierr = MPI_Allreduce(&fin,f,1,MPIU_REAL,MPIU_SUM,comm);CHKERRMPI(ierr);
256aad13602SShrirang Abhyankar   ierr = VecAssemblyBegin(G);CHKERRQ(ierr);
257aad13602SShrirang Abhyankar   ierr = VecAssemblyEnd(G);CHKERRQ(ierr);
258aad13602SShrirang Abhyankar   PetscFunctionReturn(0);
259aad13602SShrirang Abhyankar }
260aad13602SShrirang Abhyankar 
261aad13602SShrirang Abhyankar PetscErrorCode FormHessian(Tao tao, Vec x,Mat H, Mat Hpre, void *ctx)
262aad13602SShrirang Abhyankar {
263aad13602SShrirang Abhyankar   AppCtx            *user=(AppCtx*)ctx;
264aad13602SShrirang Abhyankar   Vec               DE,DI;
265aad13602SShrirang Abhyankar   const PetscScalar *de, *di;
266aad13602SShrirang Abhyankar   PetscInt          zero=0,one=1;
267aad13602SShrirang Abhyankar   PetscScalar       two=2.0;
268aad13602SShrirang Abhyankar   PetscScalar       val=0.0;
269aad13602SShrirang Abhyankar   Vec               Deseq,Diseq=user->Xseq;
270aad13602SShrirang Abhyankar   VecScatter        Descat,Discat=user->scat;
271aad13602SShrirang Abhyankar   PetscMPIInt       rank;
272aad13602SShrirang Abhyankar   MPI_Comm          comm;
273aad13602SShrirang Abhyankar   PetscErrorCode    ierr;
274aad13602SShrirang Abhyankar 
275aad13602SShrirang Abhyankar   PetscFunctionBegin;
276aad13602SShrirang Abhyankar   ierr = TaoGetDualVariables(tao,&DE,&DI);CHKERRQ(ierr);
277aad13602SShrirang Abhyankar 
278aad13602SShrirang Abhyankar   ierr = PetscObjectGetComm((PetscObject)tao,&comm);CHKERRQ(ierr);
279ffc4695bSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr);
280aad13602SShrirang Abhyankar 
2818e58fa1dSresundermann   if (!user->noeqflag){
282aad13602SShrirang Abhyankar    ierr = VecScatterCreateToZero(DE,&Descat,&Deseq);CHKERRQ(ierr);
283aad13602SShrirang Abhyankar    ierr = VecScatterBegin(Descat,DE,Deseq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
284aad13602SShrirang Abhyankar    ierr = VecScatterEnd(Descat,DE,Deseq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
2858e58fa1dSresundermann   }
2868e58fa1dSresundermann 
2878e58fa1dSresundermann   ierr = VecScatterBegin(Discat,DI,Diseq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
288aad13602SShrirang Abhyankar   ierr = VecScatterEnd(Discat,DI,Diseq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
289aad13602SShrirang Abhyankar 
290aad13602SShrirang Abhyankar   if (!rank){
2918e58fa1dSresundermann     if (!user->noeqflag){
292aad13602SShrirang Abhyankar       ierr = VecGetArrayRead(Deseq,&de);CHKERRQ(ierr);  /* places equality constraint dual into array */
2938e58fa1dSresundermann     }
2948e58fa1dSresundermann 
295aad13602SShrirang Abhyankar     ierr = VecGetArrayRead(Diseq,&di);CHKERRQ(ierr);  /* places inequality constraint dual into array */
2968e58fa1dSresundermann     if (!user->noeqflag){
297aad13602SShrirang Abhyankar       val = 2.0 * (1 + de[0] + di[0] - di[1]);
298aad13602SShrirang Abhyankar       ierr = VecRestoreArrayRead(Deseq,&de);CHKERRQ(ierr);
2998e58fa1dSresundermann     }else{
3008e58fa1dSresundermann       val = 2.0 * (1 + di[0] - di[1]);
3018e58fa1dSresundermann     }
302aad13602SShrirang Abhyankar     ierr = VecRestoreArrayRead(Diseq,&di);CHKERRQ(ierr);
303aad13602SShrirang Abhyankar     ierr = MatSetValues(H,1,&zero,1,&zero,&val,INSERT_VALUES);CHKERRQ(ierr);
304aad13602SShrirang Abhyankar     ierr = MatSetValues(H,1,&one,1,&one,&two,INSERT_VALUES);CHKERRQ(ierr);
305aad13602SShrirang Abhyankar   }
306aad13602SShrirang Abhyankar 
307aad13602SShrirang Abhyankar   ierr = MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
308aad13602SShrirang Abhyankar   ierr = MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3098e58fa1dSresundermann   if (!user->noeqflag){
310aad13602SShrirang Abhyankar     ierr = VecScatterDestroy(&Descat);CHKERRQ(ierr);
311aad13602SShrirang Abhyankar     ierr = VecDestroy(&Deseq);CHKERRQ(ierr);
3128e58fa1dSresundermann   }
313aad13602SShrirang Abhyankar   PetscFunctionReturn(0);
314aad13602SShrirang Abhyankar }
315aad13602SShrirang Abhyankar 
316aad13602SShrirang Abhyankar PetscErrorCode FormInequalityConstraints(Tao tao,Vec X,Vec CI,void *ctx)
317aad13602SShrirang Abhyankar {
318aad13602SShrirang Abhyankar   const PetscScalar *x;
319aad13602SShrirang Abhyankar   PetscScalar       ci;
320aad13602SShrirang Abhyankar   PetscErrorCode    ierr;
321aad13602SShrirang Abhyankar   MPI_Comm          comm;
322aad13602SShrirang Abhyankar   PetscMPIInt       rank;
323aad13602SShrirang Abhyankar   AppCtx            *user=(AppCtx*)ctx;
324aad13602SShrirang Abhyankar   Vec               Xseq=user->Xseq;
325aad13602SShrirang Abhyankar   VecScatter        scat=user->scat;
326aad13602SShrirang Abhyankar 
327aad13602SShrirang Abhyankar   PetscFunctionBegin;
328aad13602SShrirang Abhyankar   ierr = PetscObjectGetComm((PetscObject)tao,&comm);CHKERRQ(ierr);
329ffc4695bSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr);
330aad13602SShrirang Abhyankar 
331aad13602SShrirang Abhyankar   ierr = VecScatterBegin(scat,X,Xseq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
332aad13602SShrirang Abhyankar   ierr = VecScatterEnd(scat,X,Xseq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
333aad13602SShrirang Abhyankar 
334aad13602SShrirang Abhyankar   if (!rank) {
335aad13602SShrirang Abhyankar     ierr = VecGetArrayRead(Xseq,&x);CHKERRQ(ierr);
336aad13602SShrirang Abhyankar     ci = x[0]*x[0] - x[1];
337aad13602SShrirang Abhyankar     ierr = VecSetValue(CI,0,ci,INSERT_VALUES);CHKERRQ(ierr);
338aad13602SShrirang Abhyankar     ci = -x[0]*x[0] + x[1] + 1.0;
339aad13602SShrirang Abhyankar     ierr = VecSetValue(CI,1,ci,INSERT_VALUES);CHKERRQ(ierr);
340aad13602SShrirang Abhyankar     ierr = VecRestoreArrayRead(Xseq,&x);CHKERRQ(ierr);
341aad13602SShrirang Abhyankar   }
342aad13602SShrirang Abhyankar   ierr = VecAssemblyBegin(CI);CHKERRQ(ierr);
343aad13602SShrirang Abhyankar   ierr = VecAssemblyEnd(CI);CHKERRQ(ierr);
344aad13602SShrirang Abhyankar   PetscFunctionReturn(0);
345aad13602SShrirang Abhyankar }
346aad13602SShrirang Abhyankar 
347aad13602SShrirang Abhyankar PetscErrorCode FormEqualityConstraints(Tao tao,Vec X,Vec CE,void *ctx)
348aad13602SShrirang Abhyankar {
349aad13602SShrirang Abhyankar   const PetscScalar *x;
350aad13602SShrirang Abhyankar   PetscScalar       ce;
351aad13602SShrirang Abhyankar   PetscErrorCode    ierr;
352aad13602SShrirang Abhyankar   MPI_Comm          comm;
353aad13602SShrirang Abhyankar   PetscMPIInt       rank;
354aad13602SShrirang Abhyankar   AppCtx            *user=(AppCtx*)ctx;
355aad13602SShrirang Abhyankar   Vec               Xseq=user->Xseq;
356aad13602SShrirang Abhyankar   VecScatter        scat=user->scat;
357aad13602SShrirang Abhyankar 
358aad13602SShrirang Abhyankar   PetscFunctionBegin;
359aad13602SShrirang Abhyankar   ierr = PetscObjectGetComm((PetscObject)tao,&comm);CHKERRQ(ierr);
360ffc4695bSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr);
361aad13602SShrirang Abhyankar 
362aad13602SShrirang Abhyankar   ierr = VecScatterBegin(scat,X,Xseq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
363aad13602SShrirang Abhyankar   ierr = VecScatterEnd(scat,X,Xseq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
364aad13602SShrirang Abhyankar 
365aad13602SShrirang Abhyankar   if (!rank) {
366aad13602SShrirang Abhyankar     ierr = VecGetArrayRead(Xseq,&x);CHKERRQ(ierr);
367aad13602SShrirang Abhyankar     ce = x[0]*x[0] + x[1] - 2.0;
368aad13602SShrirang Abhyankar     ierr = VecSetValue(CE,0,ce,INSERT_VALUES);CHKERRQ(ierr);
369aad13602SShrirang Abhyankar     ierr = VecRestoreArrayRead(Xseq,&x);CHKERRQ(ierr);
370aad13602SShrirang Abhyankar   }
371aad13602SShrirang Abhyankar   ierr = VecAssemblyBegin(CE);CHKERRQ(ierr);
372aad13602SShrirang Abhyankar   ierr = VecAssemblyEnd(CE);CHKERRQ(ierr);
373aad13602SShrirang Abhyankar   PetscFunctionReturn(0);
374aad13602SShrirang Abhyankar }
375aad13602SShrirang Abhyankar 
376aad13602SShrirang Abhyankar PetscErrorCode FormInequalityJacobian(Tao tao, Vec X, Mat JI, Mat JIpre,  void *ctx)
377aad13602SShrirang Abhyankar {
378aad13602SShrirang Abhyankar   AppCtx            *user=(AppCtx*)ctx;
379aad13602SShrirang Abhyankar   PetscInt          cols[2],min,max,i;
380aad13602SShrirang Abhyankar   PetscScalar       vals[2];
381aad13602SShrirang Abhyankar   const PetscScalar *x;
382aad13602SShrirang Abhyankar   PetscErrorCode    ierr;
383aad13602SShrirang Abhyankar   Vec               Xseq=user->Xseq;
384aad13602SShrirang Abhyankar   VecScatter        scat=user->scat;
385aad13602SShrirang Abhyankar   MPI_Comm          comm;
386aad13602SShrirang Abhyankar   PetscMPIInt       rank;
387aad13602SShrirang Abhyankar 
388aad13602SShrirang Abhyankar   PetscFunctionBegin;
389aad13602SShrirang Abhyankar   ierr = PetscObjectGetComm((PetscObject)tao,&comm);CHKERRQ(ierr);
390ffc4695bSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr);
391aad13602SShrirang Abhyankar   ierr = VecScatterBegin(scat,X,Xseq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
392aad13602SShrirang Abhyankar   ierr = VecScatterEnd(scat,X,Xseq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
393aad13602SShrirang Abhyankar 
394aad13602SShrirang Abhyankar   ierr = VecGetArrayRead(Xseq,&x);CHKERRQ(ierr);
395aad13602SShrirang Abhyankar   ierr = MatGetOwnershipRange(JI,&min,&max);CHKERRQ(ierr);
396aad13602SShrirang Abhyankar 
397aad13602SShrirang Abhyankar   cols[0] = 0; cols[1] = 1;
398aad13602SShrirang Abhyankar   for (i=min;i<max;i++) {
399aad13602SShrirang Abhyankar     if (i==0){
400aad13602SShrirang Abhyankar       vals[0] = +2*x[0]; vals[1] = -1.0;
401aad13602SShrirang Abhyankar       ierr = MatSetValues(JI,1,&i,2,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
402aad13602SShrirang Abhyankar     }
403aad13602SShrirang Abhyankar     if (i==1) {
404aad13602SShrirang Abhyankar       vals[0] = -2*x[0]; vals[1] = +1.0;
405aad13602SShrirang Abhyankar       ierr = MatSetValues(JI,1,&i,2,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
406aad13602SShrirang Abhyankar     }
407aad13602SShrirang Abhyankar   }
408aad13602SShrirang Abhyankar   ierr = VecRestoreArrayRead(Xseq,&x);CHKERRQ(ierr);
409aad13602SShrirang Abhyankar 
410aad13602SShrirang Abhyankar   ierr = MatAssemblyBegin(JI,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
411aad13602SShrirang Abhyankar   ierr = MatAssemblyEnd(JI,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
412aad13602SShrirang Abhyankar   PetscFunctionReturn(0);
413aad13602SShrirang Abhyankar }
414aad13602SShrirang Abhyankar 
415aad13602SShrirang Abhyankar PetscErrorCode FormEqualityJacobian(Tao tao,Vec X,Mat JE,Mat JEpre,void *ctx)
416aad13602SShrirang Abhyankar {
417aad13602SShrirang Abhyankar   PetscInt          rows[2];
418aad13602SShrirang Abhyankar   PetscScalar       vals[2];
419aad13602SShrirang Abhyankar   const PetscScalar *x;
420aad13602SShrirang Abhyankar   PetscMPIInt       rank;
421aad13602SShrirang Abhyankar   MPI_Comm          comm;
422aad13602SShrirang Abhyankar   PetscErrorCode    ierr;
423aad13602SShrirang Abhyankar 
424aad13602SShrirang Abhyankar   PetscFunctionBegin;
425aad13602SShrirang Abhyankar   ierr = PetscObjectGetComm((PetscObject)tao,&comm);CHKERRQ(ierr);
426ffc4695bSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr);
427aad13602SShrirang Abhyankar 
428aad13602SShrirang Abhyankar   if (!rank) {
429aad13602SShrirang Abhyankar     ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
430aad13602SShrirang Abhyankar     rows[0] = 0;       rows[1] = 1;
431aad13602SShrirang Abhyankar     vals[0] = 2*x[0];  vals[1] = 1.0;
432aad13602SShrirang Abhyankar     ierr = MatSetValues(JE,1,rows,2,rows,vals,INSERT_VALUES);CHKERRQ(ierr);
433aad13602SShrirang Abhyankar     ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
434aad13602SShrirang Abhyankar   }
435aad13602SShrirang Abhyankar   ierr = MatAssemblyBegin(JE,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
436aad13602SShrirang Abhyankar   ierr = MatAssemblyEnd(JE,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
437aad13602SShrirang Abhyankar   PetscFunctionReturn(0);
438aad13602SShrirang Abhyankar }
439aad13602SShrirang Abhyankar 
440aad13602SShrirang Abhyankar 
441aad13602SShrirang Abhyankar /*TEST
442aad13602SShrirang Abhyankar 
443aad13602SShrirang Abhyankar    build:
44412d688e0SRylee Sundermann       requires: !complex !define(PETSC_USE_CXX) mumps
445aad13602SShrirang Abhyankar 
446aad13602SShrirang Abhyankar    test:
447aad13602SShrirang Abhyankar       args: -tao_converged_reason
448aad13602SShrirang Abhyankar 
449aad13602SShrirang Abhyankar    test:
450aad13602SShrirang Abhyankar       suffix: 2
451aad13602SShrirang Abhyankar       nsize: 2
452aad13602SShrirang Abhyankar       args: -tao_converged_reason
453aad13602SShrirang Abhyankar 
4548e58fa1dSresundermann    test:
4558e58fa1dSresundermann       suffix: 3
4568e58fa1dSresundermann       args: -tao_converged_reason -no_eq
4578e58fa1dSresundermann 
4588e58fa1dSresundermann    test:
4598e58fa1dSresundermann       suffix: 4
4608e58fa1dSresundermann       nsize: 2
4618e58fa1dSresundermann       args: -tao_converged_reason -no_eq
4628e58fa1dSresundermann 
463*661095bbSAlp Dener    test:
464*661095bbSAlp Dener       suffix: 5
465*661095bbSAlp Dener       args: -tao_cmonitor -tao_type almm
466*661095bbSAlp Dener 
467*661095bbSAlp Dener    test:
468*661095bbSAlp Dener       suffix: 6
469*661095bbSAlp Dener       args: -tao_cmonitor -tao_type almm -tao_almm_type phr
470*661095bbSAlp Dener 
471*661095bbSAlp Dener    test:
472*661095bbSAlp Dener       suffix: 7
473*661095bbSAlp Dener       nsize: 2
474*661095bbSAlp Dener       args: -tao_cmonitor -tao_type almm
475*661095bbSAlp Dener 
476*661095bbSAlp Dener    test:
477*661095bbSAlp Dener       suffix: 8
478*661095bbSAlp Dener       nsize: 2
479*661095bbSAlp Dener       requires: cuda
480*661095bbSAlp Dener       args: -tao_cmonitor -tao_type almm -vec_type cuda -mat_type aijcusparse
481*661095bbSAlp Dener 
482*661095bbSAlp Dener    test:
483*661095bbSAlp Dener       suffix: 9
484*661095bbSAlp Dener       nsize: 2
485*661095bbSAlp Dener       args: -tao_cmonitor -tao_type almm -no_eq
486*661095bbSAlp Dener 
487*661095bbSAlp Dener    test:
488*661095bbSAlp Dener       suffix: 10
489*661095bbSAlp Dener       nsize: 2
490*661095bbSAlp Dener       args: -tao_cmonitor -tao_type almm -tao_almm_type phr -no_eq
491*661095bbSAlp Dener 
492aad13602SShrirang Abhyankar TEST*/
493