xref: /petsc/src/tao/constrained/tutorials/ex1.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
1aad13602SShrirang Abhyankar /* Program usage: mpiexec -n 2 ./ex1 [-help] [all TAO options] */
2aad13602SShrirang Abhyankar 
3aad13602SShrirang Abhyankar /* ----------------------------------------------------------------------
4628da978Sresundermann min f(x) = (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
8628da978Sresundermann -->
9628da978Sresundermann       g(x)  = 0
10628da978Sresundermann       h(x) >= 0
11628da978Sresundermann       -1 <= x0, x1 <= 2
12628da978Sresundermann where
13628da978Sresundermann       g(x) = x0^2 + x1 - 2
14628da978Sresundermann       h(x) = [x0^2 - x1
15628da978Sresundermann               1 -(x0^2 - x1)]
16aad13602SShrirang Abhyankar ---------------------------------------------------------------------- */
17aad13602SShrirang Abhyankar 
18aad13602SShrirang Abhyankar #include <petsctao.h>
19aad13602SShrirang Abhyankar 
20aad13602SShrirang Abhyankar static char help[] = "Solves constrained optimiztion problem using pdipm.\n\
21aad13602SShrirang Abhyankar Input parameters include:\n\
22aad13602SShrirang Abhyankar   -tao_type pdipm    : sets Tao solver\n\
238e58fa1dSresundermann   -no_eq             : removes the equaility constraints from the problem\n\
24f560b561SHong Zhang   -init_view         : view initial object setup\n\
25aad13602SShrirang Abhyankar   -snes_fd           : snes with finite difference Jacobian (needed for pdipm)\n\
26628da978Sresundermann   -snes_compare_explicit : compare user Jacobian with finite difference Jacobian \n\
27aad13602SShrirang Abhyankar   -tao_cmonitor      : convergence monitor with constraint norm \n\
28a5b23f4aSJose E. Roman   -tao_view_solution : view exact solution at each iteration\n\
29f560b561SHong Zhang   Note: external package MUMPS is required to run pdipm in parallel. This is designed for a maximum of 2 processors, the code will error if size > 2.\n";
30aad13602SShrirang Abhyankar 
31aad13602SShrirang Abhyankar /*
32628da978Sresundermann    User-defined application context - contains data needed by the application
33aad13602SShrirang Abhyankar */
34aad13602SShrirang Abhyankar typedef struct {
35aad13602SShrirang Abhyankar   PetscInt   n;  /* Global length of x */
36aad13602SShrirang Abhyankar   PetscInt   ne; /* Global number of equality constraints */
37aad13602SShrirang Abhyankar   PetscInt   ni; /* Global number of inequality constraints */
38f560b561SHong Zhang   PetscBool  noeqflag, initview;
39aad13602SShrirang Abhyankar   Vec        x, xl, xu;
40aad13602SShrirang Abhyankar   Vec        ce, ci, bl, bu, Xseq;
41aad13602SShrirang Abhyankar   Mat        Ae, Ai, H;
42aad13602SShrirang Abhyankar   VecScatter scat;
43aad13602SShrirang Abhyankar } AppCtx;
44aad13602SShrirang Abhyankar 
45aad13602SShrirang Abhyankar /* -------- User-defined Routines --------- */
46aad13602SShrirang Abhyankar PetscErrorCode InitializeProblem(AppCtx *);
47aad13602SShrirang Abhyankar PetscErrorCode DestroyProblem(AppCtx *);
48aad13602SShrirang Abhyankar PetscErrorCode FormFunctionGradient(Tao, Vec, PetscReal *, Vec, void *);
49aad13602SShrirang Abhyankar PetscErrorCode FormHessian(Tao, Vec, Mat, Mat, void *);
50aad13602SShrirang Abhyankar PetscErrorCode FormInequalityConstraints(Tao, Vec, Vec, void *);
51aad13602SShrirang Abhyankar PetscErrorCode FormEqualityConstraints(Tao, Vec, Vec, void *);
52aad13602SShrirang Abhyankar PetscErrorCode FormInequalityJacobian(Tao, Vec, Mat, Mat, void *);
53aad13602SShrirang Abhyankar PetscErrorCode FormEqualityJacobian(Tao, Vec, Mat, Mat, void *);
54aad13602SShrirang Abhyankar 
55*d71ae5a4SJacob Faibussowitsch PetscErrorCode main(int argc, char **argv)
56*d71ae5a4SJacob Faibussowitsch {
57aad13602SShrirang Abhyankar   Tao         tao;
58aad13602SShrirang Abhyankar   KSP         ksp;
59aad13602SShrirang Abhyankar   PC          pc;
60aad13602SShrirang Abhyankar   AppCtx      user; /* application context */
61f560b561SHong Zhang   Vec         x, G, CI, CE;
62f560b561SHong Zhang   PetscMPIInt size;
63661095bbSAlp Dener   TaoType     type;
64f560b561SHong Zhang   PetscReal   f;
65661095bbSAlp Dener   PetscBool   pdipm;
66aad13602SShrirang Abhyankar 
67327415f7SBarry Smith   PetscFunctionBeginUser;
689566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
699566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
703c859ba3SBarry Smith   PetscCheck(size <= 2, PETSC_COMM_SELF, PETSC_ERR_WRONG_MPI_SIZE, "More than 2 processors detected. Example written to use max of 2 processors.");
71aad13602SShrirang Abhyankar 
729566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "---- Constrained Problem -----\n"));
739566063dSJacob Faibussowitsch   PetscCall(InitializeProblem(&user)); /* sets up problem, function below */
749566063dSJacob Faibussowitsch   PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao));
759566063dSJacob Faibussowitsch   PetscCall(TaoSetType(tao, TAOPDIPM));
769566063dSJacob Faibussowitsch   PetscCall(TaoSetSolution(tao, user.x));                 /* gets solution vector from problem */
779566063dSJacob Faibussowitsch   PetscCall(TaoSetVariableBounds(tao, user.xl, user.xu)); /* sets lower upper bounds from given solution */
789566063dSJacob Faibussowitsch   PetscCall(TaoSetObjectiveAndGradient(tao, NULL, FormFunctionGradient, (void *)&user));
79aad13602SShrirang Abhyankar 
8048a46eb9SPierre Jolivet   if (!user.noeqflag) PetscCall(TaoSetEqualityConstraintsRoutine(tao, user.ce, FormEqualityConstraints, (void *)&user));
819566063dSJacob Faibussowitsch   PetscCall(TaoSetInequalityConstraintsRoutine(tao, user.ci, FormInequalityConstraints, (void *)&user));
829371c9d4SSatish Balay   if (!user.noeqflag) { PetscCall(TaoSetJacobianEqualityRoutine(tao, user.Ae, user.Ae, FormEqualityJacobian, (void *)&user)); /* equality jacobian */ }
839566063dSJacob Faibussowitsch   PetscCall(TaoSetJacobianInequalityRoutine(tao, user.Ai, user.Ai, FormInequalityJacobian, (void *)&user)); /* inequality jacobian */
849566063dSJacob Faibussowitsch   PetscCall(TaoSetTolerances(tao, 1.e-6, 1.e-6, 1.e-6));
859566063dSJacob Faibussowitsch   PetscCall(TaoSetConstraintTolerances(tao, 1.e-6, 1.e-6));
86aad13602SShrirang Abhyankar 
879566063dSJacob Faibussowitsch   PetscCall(TaoGetKSP(tao, &ksp));
889566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(ksp, &pc));
899566063dSJacob Faibussowitsch   PetscCall(PCSetType(pc, PCCHOLESKY));
90aad13602SShrirang Abhyankar   /*
91aad13602SShrirang Abhyankar       This algorithm produces matrices with zeros along the diagonal therefore we use
9212d688e0SRylee Sundermann     MUMPS which provides solver for indefinite matrices
93aad13602SShrirang Abhyankar   */
94f560b561SHong Zhang #if defined(PETSC_HAVE_MUMPS)
959566063dSJacob Faibussowitsch   PetscCall(PCFactorSetMatSolverType(pc, MATSOLVERMUMPS)); /* requires mumps to solve pdipm */
96f560b561SHong Zhang #else
973c859ba3SBarry Smith   PetscCheck(size == 1, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Requires an external package that supports parallel PCCHOLESKY, e.g., MUMPS.");
98f560b561SHong Zhang #endif
999566063dSJacob Faibussowitsch   PetscCall(KSPSetType(ksp, KSPPREONLY));
1009566063dSJacob Faibussowitsch   PetscCall(KSPSetFromOptions(ksp));
101aad13602SShrirang Abhyankar 
1029566063dSJacob Faibussowitsch   PetscCall(TaoSetFromOptions(tao));
1039566063dSJacob Faibussowitsch   PetscCall(TaoGetType(tao, &type));
1049566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)tao, TAOPDIPM, &pdipm));
105661095bbSAlp Dener   if (pdipm) {
1069566063dSJacob Faibussowitsch     PetscCall(TaoSetHessian(tao, user.H, user.H, FormHessian, (void *)&user));
107f560b561SHong Zhang     if (user.initview) {
1089566063dSJacob Faibussowitsch       PetscCall(TaoSetUp(tao));
1099566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(user.x, &G));
1109566063dSJacob Faibussowitsch       PetscCall(FormFunctionGradient(tao, user.x, &f, G, (void *)&user));
1119566063dSJacob Faibussowitsch       PetscCall(FormHessian(tao, user.x, user.H, user.H, (void *)&user));
1129566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushTab(PETSC_VIEWER_STDOUT_WORLD));
11363a3b9bcSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nInitial point X:\n"));
1149566063dSJacob Faibussowitsch       PetscCall(VecView(user.x, PETSC_VIEWER_STDOUT_WORLD));
11563a3b9bcSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nInitial objective f(x) = %g\n", (double)f));
11663a3b9bcSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nInitial gradient and Hessian:\n"));
1179566063dSJacob Faibussowitsch       PetscCall(VecView(G, PETSC_VIEWER_STDOUT_WORLD));
1189566063dSJacob Faibussowitsch       PetscCall(MatView(user.H, PETSC_VIEWER_STDOUT_WORLD));
1199566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&G));
1209566063dSJacob Faibussowitsch       PetscCall(FormInequalityJacobian(tao, user.x, user.Ai, user.Ai, (void *)&user));
1219566063dSJacob Faibussowitsch       PetscCall(MatCreateVecs(user.Ai, NULL, &CI));
1229566063dSJacob Faibussowitsch       PetscCall(FormInequalityConstraints(tao, user.x, CI, (void *)&user));
12363a3b9bcSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nInitial inequality constraints and Jacobian:\n"));
1249566063dSJacob Faibussowitsch       PetscCall(VecView(CI, PETSC_VIEWER_STDOUT_WORLD));
1259566063dSJacob Faibussowitsch       PetscCall(MatView(user.Ai, PETSC_VIEWER_STDOUT_WORLD));
1269566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&CI));
127f560b561SHong Zhang       if (!user.noeqflag) {
1289566063dSJacob Faibussowitsch         PetscCall(FormEqualityJacobian(tao, user.x, user.Ae, user.Ae, (void *)&user));
1299566063dSJacob Faibussowitsch         PetscCall(MatCreateVecs(user.Ae, NULL, &CE));
1309566063dSJacob Faibussowitsch         PetscCall(FormEqualityConstraints(tao, user.x, CE, (void *)&user));
13163a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nInitial equality constraints and Jacobian:\n"));
1329566063dSJacob Faibussowitsch         PetscCall(VecView(CE, PETSC_VIEWER_STDOUT_WORLD));
1339566063dSJacob Faibussowitsch         PetscCall(MatView(user.Ae, PETSC_VIEWER_STDOUT_WORLD));
1349566063dSJacob Faibussowitsch         PetscCall(VecDestroy(&CE));
135f560b561SHong Zhang       }
1369566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n"));
1379566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopTab(PETSC_VIEWER_STDOUT_WORLD));
138f560b561SHong Zhang     }
139661095bbSAlp Dener   }
140aad13602SShrirang Abhyankar 
1419566063dSJacob Faibussowitsch   PetscCall(TaoSolve(tao));
1429566063dSJacob Faibussowitsch   PetscCall(TaoGetSolution(tao, &x));
1439566063dSJacob Faibussowitsch   PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));
144aad13602SShrirang Abhyankar 
145aad13602SShrirang Abhyankar   /* Free objects */
1469566063dSJacob Faibussowitsch   PetscCall(DestroyProblem(&user));
1479566063dSJacob Faibussowitsch   PetscCall(TaoDestroy(&tao));
1489566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
149b122ec5aSJacob Faibussowitsch   return 0;
150aad13602SShrirang Abhyankar }
151aad13602SShrirang Abhyankar 
152*d71ae5a4SJacob Faibussowitsch PetscErrorCode InitializeProblem(AppCtx *user)
153*d71ae5a4SJacob Faibussowitsch {
154aad13602SShrirang Abhyankar   PetscMPIInt size;
155aad13602SShrirang Abhyankar   PetscMPIInt rank;
156aad13602SShrirang Abhyankar   PetscInt    nloc, neloc, niloc;
157aad13602SShrirang Abhyankar 
158aad13602SShrirang Abhyankar   PetscFunctionBegin;
1599566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
1609566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
1618e58fa1dSresundermann   user->noeqflag = PETSC_FALSE;
1629566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-no_eq", &user->noeqflag, NULL));
163f560b561SHong Zhang   user->initview = PETSC_FALSE;
1649566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-init_view", &user->initview, NULL));
165628da978Sresundermann 
1668e58fa1dSresundermann   if (!user->noeqflag) {
167f560b561SHong Zhang     /* Tell user the correct solution, not an error checking */
1689566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Solution should be f(1,1)=-2\n"));
1698e58fa1dSresundermann   }
170aad13602SShrirang Abhyankar 
1712d4ee042Sprj-   /* create vector x and set initial values */
172aad13602SShrirang Abhyankar   user->n = 2; /* global length */
173f560b561SHong Zhang   nloc    = (size == 1) ? user->n : 1;
1749566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &user->x));
1759566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(user->x, nloc, user->n));
1769566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(user->x));
1779566063dSJacob Faibussowitsch   PetscCall(VecSet(user->x, 0));
178aad13602SShrirang Abhyankar 
179aad13602SShrirang Abhyankar   /* create and set lower and upper bound vectors */
1809566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user->x, &user->xl));
1819566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user->x, &user->xu));
1829566063dSJacob Faibussowitsch   PetscCall(VecSet(user->xl, -1.0));
1839566063dSJacob Faibussowitsch   PetscCall(VecSet(user->xu, 2.0));
184aad13602SShrirang Abhyankar 
185aad13602SShrirang Abhyankar   /* create scater to zero */
1869566063dSJacob Faibussowitsch   PetscCall(VecScatterCreateToZero(user->x, &user->scat, &user->Xseq));
187aad13602SShrirang Abhyankar 
188aad13602SShrirang Abhyankar   user->ne = 1;
189628da978Sresundermann   user->ni = 2;
190aad13602SShrirang Abhyankar   neloc    = (rank == 0) ? user->ne : 0;
191f560b561SHong Zhang   niloc    = (size == 1) ? user->ni : 1;
192aad13602SShrirang Abhyankar 
1938e58fa1dSresundermann   if (!user->noeqflag) {
1949566063dSJacob Faibussowitsch     PetscCall(VecCreate(PETSC_COMM_WORLD, &user->ce)); /* a 1x1 vec for equality constraints */
1959566063dSJacob Faibussowitsch     PetscCall(VecSetSizes(user->ce, neloc, user->ne));
1969566063dSJacob Faibussowitsch     PetscCall(VecSetFromOptions(user->ce));
1979566063dSJacob Faibussowitsch     PetscCall(VecSetUp(user->ce));
1988e58fa1dSresundermann   }
199628da978Sresundermann 
2009566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &user->ci)); /* a 2x1 vec for inequality constraints */
2019566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(user->ci, niloc, user->ni));
2029566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(user->ci));
2039566063dSJacob Faibussowitsch   PetscCall(VecSetUp(user->ci));
204aad13602SShrirang Abhyankar 
205a5b23f4aSJose E. Roman   /* nexn & nixn matricies for equally and inequalty constraints */
2068e58fa1dSresundermann   if (!user->noeqflag) {
2079566063dSJacob Faibussowitsch     PetscCall(MatCreate(PETSC_COMM_WORLD, &user->Ae));
2089566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(user->Ae, neloc, nloc, user->ne, user->n));
2099566063dSJacob Faibussowitsch     PetscCall(MatSetFromOptions(user->Ae));
2109566063dSJacob Faibussowitsch     PetscCall(MatSetUp(user->Ae));
2118e58fa1dSresundermann   }
2128e58fa1dSresundermann 
2139566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &user->Ai));
2149566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(user->Ai, niloc, nloc, user->ni, user->n));
2159566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(user->Ai));
2169566063dSJacob Faibussowitsch   PetscCall(MatSetUp(user->Ai));
217628da978Sresundermann 
2189566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &user->H));
2199566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(user->H, nloc, nloc, user->n, user->n));
2209566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(user->H));
2219566063dSJacob Faibussowitsch   PetscCall(MatSetUp(user->H));
222aad13602SShrirang Abhyankar   PetscFunctionReturn(0);
223aad13602SShrirang Abhyankar }
224aad13602SShrirang Abhyankar 
225*d71ae5a4SJacob Faibussowitsch PetscErrorCode DestroyProblem(AppCtx *user)
226*d71ae5a4SJacob Faibussowitsch {
227aad13602SShrirang Abhyankar   PetscFunctionBegin;
22848a46eb9SPierre Jolivet   if (!user->noeqflag) PetscCall(MatDestroy(&user->Ae));
2299566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->Ai));
2309566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user->H));
231aad13602SShrirang Abhyankar 
2329566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->x));
23348a46eb9SPierre Jolivet   if (!user->noeqflag) PetscCall(VecDestroy(&user->ce));
2349566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->ci));
2359566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->xl));
2369566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->xu));
2379566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user->Xseq));
2389566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&user->scat));
239aad13602SShrirang Abhyankar   PetscFunctionReturn(0);
240aad13602SShrirang Abhyankar }
241aad13602SShrirang Abhyankar 
242628da978Sresundermann /* Evaluate
243628da978Sresundermann    f(x) = (x0 - 2)^2 + (x1 - 2)^2 - 2*(x0 + x1)
244628da978Sresundermann    G = grad f = [2*(x0 - 2) - 2;
245628da978Sresundermann                  2*(x1 - 2) - 2]
246aad13602SShrirang Abhyankar */
247*d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *f, Vec G, void *ctx)
248*d71ae5a4SJacob Faibussowitsch {
249aad13602SShrirang Abhyankar   PetscScalar        g;
250aad13602SShrirang Abhyankar   const PetscScalar *x;
251aad13602SShrirang Abhyankar   MPI_Comm           comm;
252aad13602SShrirang Abhyankar   PetscMPIInt        rank;
253aad13602SShrirang Abhyankar   PetscReal          fin;
254aad13602SShrirang Abhyankar   AppCtx            *user = (AppCtx *)ctx;
255aad13602SShrirang Abhyankar   Vec                Xseq = user->Xseq;
256aad13602SShrirang Abhyankar   VecScatter         scat = user->scat;
257aad13602SShrirang Abhyankar 
258aad13602SShrirang Abhyankar   PetscFunctionBegin;
2599566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)tao, &comm));
2609566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
261aad13602SShrirang Abhyankar 
2629566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(scat, X, Xseq, INSERT_VALUES, SCATTER_FORWARD));
2639566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(scat, X, Xseq, INSERT_VALUES, SCATTER_FORWARD));
264aad13602SShrirang Abhyankar 
265aad13602SShrirang Abhyankar   fin = 0.0;
266dd400576SPatrick Sanan   if (rank == 0) {
2679566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(Xseq, &x));
268aad13602SShrirang Abhyankar     fin = (x[0] - 2.0) * (x[0] - 2.0) + (x[1] - 2.0) * (x[1] - 2.0) - 2.0 * (x[0] + x[1]);
269aad13602SShrirang Abhyankar     g   = 2.0 * (x[0] - 2.0) - 2.0;
2709566063dSJacob Faibussowitsch     PetscCall(VecSetValue(G, 0, g, INSERT_VALUES));
271aad13602SShrirang Abhyankar     g = 2.0 * (x[1] - 2.0) - 2.0;
2729566063dSJacob Faibussowitsch     PetscCall(VecSetValue(G, 1, g, INSERT_VALUES));
2739566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(Xseq, &x));
274aad13602SShrirang Abhyankar   }
2759566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allreduce(&fin, f, 1, MPIU_REAL, MPIU_SUM, comm));
2769566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(G));
2779566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(G));
278aad13602SShrirang Abhyankar   PetscFunctionReturn(0);
279aad13602SShrirang Abhyankar }
280aad13602SShrirang Abhyankar 
281628da978Sresundermann /* Evaluate
282628da978Sresundermann    H = fxx + grad (grad g^T*DI) - grad (grad h^T*DE)]
283628da978Sresundermann      = [ 2*(1+de[0]-di[0]+di[1]), 0;
284628da978Sresundermann                    0,             2]
285628da978Sresundermann */
286*d71ae5a4SJacob Faibussowitsch PetscErrorCode FormHessian(Tao tao, Vec x, Mat H, Mat Hpre, void *ctx)
287*d71ae5a4SJacob Faibussowitsch {
288aad13602SShrirang Abhyankar   AppCtx            *user = (AppCtx *)ctx;
289aad13602SShrirang Abhyankar   Vec                DE, DI;
290aad13602SShrirang Abhyankar   const PetscScalar *de, *di;
291aad13602SShrirang Abhyankar   PetscInt           zero = 0, one = 1;
292aad13602SShrirang Abhyankar   PetscScalar        two = 2.0;
293aad13602SShrirang Abhyankar   PetscScalar        val = 0.0;
294f560b561SHong Zhang   Vec                Deseq, Diseq;
295f560b561SHong Zhang   VecScatter         Descat, Discat;
296aad13602SShrirang Abhyankar   PetscMPIInt        rank;
297aad13602SShrirang Abhyankar   MPI_Comm           comm;
298aad13602SShrirang Abhyankar 
299aad13602SShrirang Abhyankar   PetscFunctionBegin;
3009566063dSJacob Faibussowitsch   PetscCall(TaoGetDualVariables(tao, &DE, &DI));
301aad13602SShrirang Abhyankar 
3029566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)tao, &comm));
3039566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
304aad13602SShrirang Abhyankar 
3058e58fa1dSresundermann   if (!user->noeqflag) {
3069566063dSJacob Faibussowitsch     PetscCall(VecScatterCreateToZero(DE, &Descat, &Deseq));
3079566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(Descat, DE, Deseq, INSERT_VALUES, SCATTER_FORWARD));
3089566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(Descat, DE, Deseq, INSERT_VALUES, SCATTER_FORWARD));
3098e58fa1dSresundermann   }
3109566063dSJacob Faibussowitsch   PetscCall(VecScatterCreateToZero(DI, &Discat, &Diseq));
3119566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(Discat, DI, Diseq, INSERT_VALUES, SCATTER_FORWARD));
3129566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(Discat, DI, Diseq, INSERT_VALUES, SCATTER_FORWARD));
313aad13602SShrirang Abhyankar 
314dd400576SPatrick Sanan   if (rank == 0) {
3159371c9d4SSatish Balay     if (!user->noeqflag) { PetscCall(VecGetArrayRead(Deseq, &de)); /* places equality constraint dual into array */ }
3169566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(Diseq, &di)); /* places inequality constraint dual into array */
317628da978Sresundermann 
3188e58fa1dSresundermann     if (!user->noeqflag) {
319628da978Sresundermann       val = 2.0 * (1 + de[0] - di[0] + di[1]);
3209566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(Deseq, &de));
3219566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(Diseq, &di));
3228e58fa1dSresundermann     } else {
323628da978Sresundermann       val = 2.0 * (1 - di[0] + di[1]);
3248e58fa1dSresundermann     }
3259566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(Diseq, &di));
3269566063dSJacob Faibussowitsch     PetscCall(MatSetValues(H, 1, &zero, 1, &zero, &val, INSERT_VALUES));
3279566063dSJacob Faibussowitsch     PetscCall(MatSetValues(H, 1, &one, 1, &one, &two, INSERT_VALUES));
328aad13602SShrirang Abhyankar   }
3299566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(H, MAT_FINAL_ASSEMBLY));
3309566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(H, MAT_FINAL_ASSEMBLY));
3318e58fa1dSresundermann   if (!user->noeqflag) {
3329566063dSJacob Faibussowitsch     PetscCall(VecScatterDestroy(&Descat));
3339566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&Deseq));
3348e58fa1dSresundermann   }
3359566063dSJacob Faibussowitsch   PetscCall(VecScatterDestroy(&Discat));
3369566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Diseq));
337aad13602SShrirang Abhyankar   PetscFunctionReturn(0);
338aad13602SShrirang Abhyankar }
339aad13602SShrirang Abhyankar 
340628da978Sresundermann /* Evaluate
341628da978Sresundermann    h = [ x0^2 - x1;
342628da978Sresundermann          1 -(x0^2 - x1)]
343628da978Sresundermann */
344*d71ae5a4SJacob Faibussowitsch PetscErrorCode FormInequalityConstraints(Tao tao, Vec X, Vec CI, void *ctx)
345*d71ae5a4SJacob Faibussowitsch {
346aad13602SShrirang Abhyankar   const PetscScalar *x;
347aad13602SShrirang Abhyankar   PetscScalar        ci;
348aad13602SShrirang Abhyankar   MPI_Comm           comm;
349aad13602SShrirang Abhyankar   PetscMPIInt        rank;
350aad13602SShrirang Abhyankar   AppCtx            *user = (AppCtx *)ctx;
351aad13602SShrirang Abhyankar   Vec                Xseq = user->Xseq;
352aad13602SShrirang Abhyankar   VecScatter         scat = user->scat;
353aad13602SShrirang Abhyankar 
354aad13602SShrirang Abhyankar   PetscFunctionBegin;
3559566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)tao, &comm));
3569566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
357aad13602SShrirang Abhyankar 
3589566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(scat, X, Xseq, INSERT_VALUES, SCATTER_FORWARD));
3599566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(scat, X, Xseq, INSERT_VALUES, SCATTER_FORWARD));
360aad13602SShrirang Abhyankar 
361dd400576SPatrick Sanan   if (rank == 0) {
3629566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(Xseq, &x));
363aad13602SShrirang Abhyankar     ci = x[0] * x[0] - x[1];
3649566063dSJacob Faibussowitsch     PetscCall(VecSetValue(CI, 0, ci, INSERT_VALUES));
365aad13602SShrirang Abhyankar     ci = -x[0] * x[0] + x[1] + 1.0;
3669566063dSJacob Faibussowitsch     PetscCall(VecSetValue(CI, 1, ci, INSERT_VALUES));
3679566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(Xseq, &x));
368aad13602SShrirang Abhyankar   }
3699566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(CI));
3709566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(CI));
371aad13602SShrirang Abhyankar   PetscFunctionReturn(0);
372aad13602SShrirang Abhyankar }
373aad13602SShrirang Abhyankar 
374628da978Sresundermann /* Evaluate
375628da978Sresundermann    g = [ x0^2 + x1 - 2]
376628da978Sresundermann */
377*d71ae5a4SJacob Faibussowitsch PetscErrorCode FormEqualityConstraints(Tao tao, Vec X, Vec CE, void *ctx)
378*d71ae5a4SJacob Faibussowitsch {
379aad13602SShrirang Abhyankar   const PetscScalar *x;
380aad13602SShrirang Abhyankar   PetscScalar        ce;
381aad13602SShrirang Abhyankar   MPI_Comm           comm;
382aad13602SShrirang Abhyankar   PetscMPIInt        rank;
383aad13602SShrirang Abhyankar   AppCtx            *user = (AppCtx *)ctx;
384aad13602SShrirang Abhyankar   Vec                Xseq = user->Xseq;
385aad13602SShrirang Abhyankar   VecScatter         scat = user->scat;
386aad13602SShrirang Abhyankar 
387aad13602SShrirang Abhyankar   PetscFunctionBegin;
3889566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)tao, &comm));
3899566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
390aad13602SShrirang Abhyankar 
3919566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(scat, X, Xseq, INSERT_VALUES, SCATTER_FORWARD));
3929566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(scat, X, Xseq, INSERT_VALUES, SCATTER_FORWARD));
393aad13602SShrirang Abhyankar 
394dd400576SPatrick Sanan   if (rank == 0) {
3959566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(Xseq, &x));
396aad13602SShrirang Abhyankar     ce = x[0] * x[0] + x[1] - 2.0;
3979566063dSJacob Faibussowitsch     PetscCall(VecSetValue(CE, 0, ce, INSERT_VALUES));
3989566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(Xseq, &x));
399aad13602SShrirang Abhyankar   }
4009566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(CE));
4019566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(CE));
402aad13602SShrirang Abhyankar   PetscFunctionReturn(0);
403aad13602SShrirang Abhyankar }
404aad13602SShrirang Abhyankar 
405628da978Sresundermann /*
406628da978Sresundermann   grad h = [  2*x0, -1;
407628da978Sresundermann              -2*x0,  1]
408628da978Sresundermann */
409*d71ae5a4SJacob Faibussowitsch PetscErrorCode FormInequalityJacobian(Tao tao, Vec X, Mat JI, Mat JIpre, void *ctx)
410*d71ae5a4SJacob Faibussowitsch {
411aad13602SShrirang Abhyankar   AppCtx            *user = (AppCtx *)ctx;
412f560b561SHong Zhang   PetscInt           zero = 0, one = 1, cols[2];
413aad13602SShrirang Abhyankar   PetscScalar        vals[2];
414aad13602SShrirang Abhyankar   const PetscScalar *x;
415aad13602SShrirang Abhyankar   Vec                Xseq = user->Xseq;
416aad13602SShrirang Abhyankar   VecScatter         scat = user->scat;
417aad13602SShrirang Abhyankar   MPI_Comm           comm;
418aad13602SShrirang Abhyankar   PetscMPIInt        rank;
419aad13602SShrirang Abhyankar 
420aad13602SShrirang Abhyankar   PetscFunctionBegin;
4219566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)tao, &comm));
4229566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
4239566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(scat, X, Xseq, INSERT_VALUES, SCATTER_FORWARD));
4249566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(scat, X, Xseq, INSERT_VALUES, SCATTER_FORWARD));
425aad13602SShrirang Abhyankar 
4269566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xseq, &x));
427c5853193SPierre Jolivet   if (rank == 0) {
4289371c9d4SSatish Balay     cols[0] = 0;
4299371c9d4SSatish Balay     cols[1] = 1;
4309371c9d4SSatish Balay     vals[0] = 2 * x[0];
4319371c9d4SSatish Balay     vals[1] = -1.0;
4329566063dSJacob Faibussowitsch     PetscCall(MatSetValues(JI, 1, &zero, 2, cols, vals, INSERT_VALUES));
4339371c9d4SSatish Balay     vals[0] = -2 * x[0];
4349371c9d4SSatish Balay     vals[1] = 1.0;
4359566063dSJacob Faibussowitsch     PetscCall(MatSetValues(JI, 1, &one, 2, cols, vals, INSERT_VALUES));
436aad13602SShrirang Abhyankar   }
4379566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xseq, &x));
4389566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(JI, MAT_FINAL_ASSEMBLY));
4399566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(JI, MAT_FINAL_ASSEMBLY));
440aad13602SShrirang Abhyankar   PetscFunctionReturn(0);
441aad13602SShrirang Abhyankar }
442aad13602SShrirang Abhyankar 
443628da978Sresundermann /*
444628da978Sresundermann   grad g = [2*x0
445628da978Sresundermann              1.0 ]
446628da978Sresundermann */
447*d71ae5a4SJacob Faibussowitsch PetscErrorCode FormEqualityJacobian(Tao tao, Vec X, Mat JE, Mat JEpre, void *ctx)
448*d71ae5a4SJacob Faibussowitsch {
449f560b561SHong Zhang   PetscInt           zero = 0, cols[2];
450aad13602SShrirang Abhyankar   PetscScalar        vals[2];
451aad13602SShrirang Abhyankar   const PetscScalar *x;
452aad13602SShrirang Abhyankar   PetscMPIInt        rank;
453aad13602SShrirang Abhyankar   MPI_Comm           comm;
454aad13602SShrirang Abhyankar 
455aad13602SShrirang Abhyankar   PetscFunctionBegin;
4569566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)tao, &comm));
4579566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
458aad13602SShrirang Abhyankar 
459dd400576SPatrick Sanan   if (rank == 0) {
4609566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(X, &x));
4619371c9d4SSatish Balay     cols[0] = 0;
4629371c9d4SSatish Balay     cols[1] = 1;
4639371c9d4SSatish Balay     vals[0] = 2 * x[0];
4649371c9d4SSatish Balay     vals[1] = 1.0;
4659566063dSJacob Faibussowitsch     PetscCall(MatSetValues(JE, 1, &zero, 2, cols, vals, INSERT_VALUES));
4669566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(X, &x));
467aad13602SShrirang Abhyankar   }
4689566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(JE, MAT_FINAL_ASSEMBLY));
4699566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(JE, MAT_FINAL_ASSEMBLY));
470aad13602SShrirang Abhyankar   PetscFunctionReturn(0);
471aad13602SShrirang Abhyankar }
472aad13602SShrirang Abhyankar 
473aad13602SShrirang Abhyankar /*TEST
474aad13602SShrirang Abhyankar 
475aad13602SShrirang Abhyankar    build:
476f560b561SHong Zhang       requires: !complex !defined(PETSC_USE_CXX)
477aad13602SShrirang Abhyankar 
478aad13602SShrirang Abhyankar    test:
47909ee8bb0SRylee Sundermann       args: -tao_converged_reason -tao_pdipm_kkt_shift_pd
480f560b561SHong Zhang       requires: mumps
481aad13602SShrirang Abhyankar 
482aad13602SShrirang Abhyankar    test:
483aad13602SShrirang Abhyankar       suffix: 2
484aad13602SShrirang Abhyankar       nsize: 2
48509ee8bb0SRylee Sundermann       args: -tao_converged_reason -tao_pdipm_kkt_shift_pd
486f560b561SHong Zhang       requires: mumps
487aad13602SShrirang Abhyankar 
4888e58fa1dSresundermann    test:
4898e58fa1dSresundermann       suffix: 3
4908e58fa1dSresundermann       args: -tao_converged_reason -no_eq
491f560b561SHong Zhang       requires: mumps
4928e58fa1dSresundermann 
4938e58fa1dSresundermann    test:
4948e58fa1dSresundermann       suffix: 4
4958e58fa1dSresundermann       nsize: 2
4968e58fa1dSresundermann       args: -tao_converged_reason -no_eq
497f560b561SHong Zhang       requires: mumps
4988e58fa1dSresundermann 
499661095bbSAlp Dener    test:
500661095bbSAlp Dener       suffix: 5
501661095bbSAlp Dener       args: -tao_cmonitor -tao_type almm
502f560b561SHong Zhang       requires: mumps
503661095bbSAlp Dener 
504661095bbSAlp Dener    test:
505661095bbSAlp Dener       suffix: 6
506661095bbSAlp Dener       args: -tao_cmonitor -tao_type almm -tao_almm_type phr
507f560b561SHong Zhang       requires: mumps
508661095bbSAlp Dener 
509661095bbSAlp Dener    test:
510661095bbSAlp Dener       suffix: 7
511661095bbSAlp Dener       nsize: 2
512f560b561SHong Zhang       requires: mumps
513661095bbSAlp Dener       args: -tao_cmonitor -tao_type almm
514661095bbSAlp Dener 
515661095bbSAlp Dener    test:
516661095bbSAlp Dener       suffix: 8
517661095bbSAlp Dener       nsize: 2
518f560b561SHong Zhang       requires: cuda mumps
519661095bbSAlp Dener       args: -tao_cmonitor -tao_type almm -vec_type cuda -mat_type aijcusparse
520661095bbSAlp Dener 
521661095bbSAlp Dener    test:
522661095bbSAlp Dener       suffix: 9
523661095bbSAlp Dener       nsize: 2
524661095bbSAlp Dener       args: -tao_cmonitor -tao_type almm -no_eq
525f560b561SHong Zhang       requires: mumps
526661095bbSAlp Dener 
527661095bbSAlp Dener    test:
528661095bbSAlp Dener       suffix: 10
529661095bbSAlp Dener       nsize: 2
530661095bbSAlp Dener       args: -tao_cmonitor -tao_type almm -tao_almm_type phr -no_eq
531f560b561SHong Zhang       requires: mumps
532661095bbSAlp Dener 
533aad13602SShrirang Abhyankar TEST*/
534