Actual source code: ex1.c

  1: /* Program usage: mpiexec -n 2 ./ex1 [-help] [all TAO options] */

  3: /* ----------------------------------------------------------------------
  4: min f(x) = (x0 - 2)^2 + (x1 - 2)^2 - 2*(x0 + x1)
  5: s.t.  x0^2 + x1 - 2 = 0
  6:       0  <= x0^2 - x1 <= 1
  7:       -1 <= x0, x1 <= 2
  8: -->
  9:       g(x)  = 0
 10:       h(x) >= 0
 11:       -1 <= x0, x1 <= 2
 12: where
 13:       g(x) = x0^2 + x1 - 2
 14:       h(x) = [x0^2 - x1
 15:               1 -(x0^2 - x1)]
 16: ---------------------------------------------------------------------- */

 18: #include <petsctao.h>

 20: static char help[] = "Solves constrained optimization problem using pdipm.\n\
 21: Input parameters include:\n\
 22:   -tao_type pdipm    : sets Tao solver\n\
 23:   -no_eq             : removes the equaility constraints from the problem\n\
 24:   -no_bound          : removes the bound constraints from the problem\n\
 25:   -init_view         : view initial object setup\n\
 26:   -snes_fd           : snes with finite difference Jacobian (needed for pdipm)\n\
 27:   -snes_compare_explicit : compare user Jacobian with finite difference Jacobian \n\
 28:   -tao_monitor_constraint_norm : convergence monitor with constraint norm \n\
 29:   -tao_monitor_solution : view exact solution at each iteration\n\
 30:   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";

 32: /*
 33:    User-defined application context - contains data needed by the application
 34: */
 35: typedef struct {
 36:   PetscInt   n;  /* Global length of x */
 37:   PetscInt   ne; /* Global number of equality constraints */
 38:   PetscInt   ni; /* Global number of inequality constraints */
 39:   PetscBool  noeqflag, noboundflag, initview;
 40:   Vec        x, xl, xu;
 41:   Vec        ce, ci, bl, bu, Xseq;
 42:   Mat        Ae, Ai, H;
 43:   VecScatter scat;
 44: } AppCtx;

 46: /* -------- User-defined Routines --------- */
 47: PetscErrorCode InitializeProblem(AppCtx *);
 48: PetscErrorCode DestroyProblem(AppCtx *);
 49: PetscErrorCode FormFunctionGradient(Tao, Vec, PetscReal *, Vec, void *);
 50: PetscErrorCode FormPDIPMHessian(Tao, Vec, Mat, Mat, void *);
 51: PetscErrorCode FormInequalityConstraints(Tao, Vec, Vec, void *);
 52: PetscErrorCode FormEqualityConstraints(Tao, Vec, Vec, void *);
 53: PetscErrorCode FormInequalityJacobian(Tao, Vec, Mat, Mat, void *);
 54: PetscErrorCode FormEqualityJacobian(Tao, Vec, Mat, Mat, void *);

 56: int main(int argc, char **argv)
 57: {
 58:   Tao         tao;
 59:   KSP         ksp;
 60:   PC          pc;
 61:   AppCtx      user; /* application context */
 62:   Vec         x, G, CI, CE;
 63:   PetscMPIInt size;
 64:   PetscReal   f;
 65:   PetscBool   pdipm;

 67:   PetscFunctionBeginUser;
 68:   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
 69:   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
 70:   PetscCheck(size <= 2, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "More than 2 processes detected. Example written to use at most 2 processes.");

 72:   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "---- Constrained Problem -----\n"));
 73:   PetscCall(InitializeProblem(&user)); /* sets up problem, function below */

 75:   PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao));
 76:   PetscCall(TaoSetType(tao, TAOALMM));
 77:   PetscCall(TaoSetSolution(tao, user.x));
 78:   if (!user.noboundflag) PetscCall(TaoSetVariableBounds(tao, user.xl, user.xu));
 79:   PetscCall(TaoSetObjectiveAndGradient(tao, NULL, FormFunctionGradient, (void *)&user));
 80:   PetscCall(TaoSetTolerances(tao, 1.e-4, 0.0, 0.0));
 81:   PetscCall(TaoSetConstraintTolerances(tao, 1.e-4, 0.0));
 82:   PetscCall(TaoSetMaximumFunctionEvaluations(tao, 1e6));
 83:   PetscCall(TaoSetFromOptions(tao));

 85:   if (!user.noeqflag) {
 86:     PetscCall(TaoSetEqualityConstraintsRoutine(tao, user.ce, FormEqualityConstraints, (void *)&user));
 87:     PetscCall(TaoSetJacobianEqualityRoutine(tao, user.Ae, user.Ae, FormEqualityJacobian, (void *)&user));
 88:   }
 89:   PetscCall(TaoSetInequalityConstraintsRoutine(tao, user.ci, FormInequalityConstraints, (void *)&user));
 90:   PetscCall(TaoSetJacobianInequalityRoutine(tao, user.Ai, user.Ai, FormInequalityJacobian, (void *)&user));

 92:   PetscCall(PetscObjectTypeCompare((PetscObject)tao, TAOPDIPM, &pdipm));
 93:   if (pdipm) {
 94:     /*
 95:       PDIPM produces an indefinite KKT matrix with some zeros along the diagonal
 96:       Inverting this indefinite matrix requires PETSc to be configured with MUMPS
 97:     */
 98:     PetscCheck(PetscDefined(HAVE_MUMPS), PetscObjectComm((PetscObject)tao), PETSC_ERR_SUP, "TAOPDIPM requires PETSc to be configured with MUMPS (--download-mumps)");
 99:     PetscCall(TaoGetKSP(tao, &ksp));
100:     PetscCall(KSPGetPC(ksp, &pc));
101:     PetscCall(PCFactorSetMatSolverType(pc, MATSOLVERMUMPS));
102:     PetscCall(KSPSetType(ksp, KSPPREONLY));
103:     PetscCall(PCSetType(pc, PCCHOLESKY));
104:     PetscCall(KSPSetFromOptions(ksp));
105:     PetscCall(TaoSetHessian(tao, user.H, user.H, FormPDIPMHessian, (void *)&user));
106:   }

108:   /* Print out an initial view of the problem */
109:   if (user.initview) {
110:     PetscCall(TaoSetUp(tao));
111:     PetscCall(VecDuplicate(user.x, &G));
112:     PetscCall(FormFunctionGradient(tao, user.x, &f, G, (void *)&user));
113:     if (pdipm) PetscCall(FormPDIPMHessian(tao, user.x, user.H, user.H, (void *)&user));
114:     PetscCall(PetscViewerASCIIPushTab(PETSC_VIEWER_STDOUT_WORLD));
115:     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nInitial point X:\n"));
116:     PetscCall(VecView(user.x, PETSC_VIEWER_STDOUT_WORLD));
117:     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nInitial objective f(x) = %g\n", (double)f));
118:     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nInitial gradient and Hessian:\n"));
119:     PetscCall(VecView(G, PETSC_VIEWER_STDOUT_WORLD));
120:     if (pdipm) PetscCall(MatView(user.H, PETSC_VIEWER_STDOUT_WORLD));
121:     PetscCall(VecDestroy(&G));
122:     PetscCall(FormInequalityJacobian(tao, user.x, user.Ai, user.Ai, (void *)&user));
123:     PetscCall(MatCreateVecs(user.Ai, NULL, &CI));
124:     PetscCall(FormInequalityConstraints(tao, user.x, CI, (void *)&user));
125:     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nInitial inequality constraints and Jacobian:\n"));
126:     PetscCall(VecView(CI, PETSC_VIEWER_STDOUT_WORLD));
127:     PetscCall(MatView(user.Ai, PETSC_VIEWER_STDOUT_WORLD));
128:     PetscCall(VecDestroy(&CI));
129:     if (!user.noeqflag) {
130:       PetscCall(FormEqualityJacobian(tao, user.x, user.Ae, user.Ae, (void *)&user));
131:       PetscCall(MatCreateVecs(user.Ae, NULL, &CE));
132:       PetscCall(FormEqualityConstraints(tao, user.x, CE, (void *)&user));
133:       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nInitial equality constraints and Jacobian:\n"));
134:       PetscCall(VecView(CE, PETSC_VIEWER_STDOUT_WORLD));
135:       PetscCall(MatView(user.Ae, PETSC_VIEWER_STDOUT_WORLD));
136:       PetscCall(VecDestroy(&CE));
137:     }
138:     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n"));
139:     PetscCall(PetscViewerASCIIPopTab(PETSC_VIEWER_STDOUT_WORLD));
140:   }

142:   PetscCall(TaoSolve(tao));
143:   PetscCall(TaoGetSolution(tao, &x));
144:   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Found solution:\n"));
145:   PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));

147:   /* Free objects */
148:   PetscCall(DestroyProblem(&user));
149:   PetscCall(TaoDestroy(&tao));
150:   PetscCall(PetscFinalize());
151:   return PETSC_SUCCESS;
152: }

154: PetscErrorCode InitializeProblem(AppCtx *user)
155: {
156:   PetscMPIInt size;
157:   PetscMPIInt rank;
158:   PetscInt    nloc, neloc, niloc;

160:   PetscFunctionBegin;
161:   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
162:   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
163:   user->noeqflag    = PETSC_FALSE;
164:   user->noboundflag = PETSC_FALSE;
165:   PetscCall(PetscOptionsGetBool(NULL, NULL, "-no_eq", &user->noeqflag, NULL));
166:   PetscCall(PetscOptionsGetBool(NULL, NULL, "-no_bound", &user->noboundflag, NULL));
167:   user->initview = PETSC_FALSE;
168:   PetscCall(PetscOptionsGetBool(NULL, NULL, "-init_view", &user->initview, NULL));

170:   /* Tell user the correct solution, not an error checking */
171:   if (!user->noeqflag) {
172:     /* With equality constraint, bound or no bound makes no different to the solution */
173:     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Expected solution: f(1, 1) = -2\n"));
174:   } else {
175:     if (user->noboundflag) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Expected solution (-no_eq, -no_bound): f(2.05655, 3.22938) = -9.05728\n"));
176:     else PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Expected solution (-no_eq): f(1.73205, 2) = -7.3923\n"));
177:   }

179:   /* create vector x and set initial values */
180:   user->n = 2; /* global length */
181:   nloc    = (size == 1) ? user->n : 1;
182:   PetscCall(VecCreate(PETSC_COMM_WORLD, &user->x));
183:   PetscCall(VecSetSizes(user->x, nloc, user->n));
184:   PetscCall(VecSetFromOptions(user->x));

186:   /* create and set lower and upper bound vectors */
187:   PetscCall(VecDuplicate(user->x, &user->xl));
188:   PetscCall(VecDuplicate(user->x, &user->xu));
189:   PetscCall(VecSet(user->xl, -1.0));
190:   PetscCall(VecSet(user->xu, 2.0));

192:   /* create scater to zero */
193:   PetscCall(VecScatterCreateToZero(user->x, &user->scat, &user->Xseq));

195:   user->ne = 1;
196:   user->ni = 2;
197:   neloc    = (rank == 0) ? user->ne : 0;
198:   niloc    = (size == 1) ? user->ni : 1;

200:   if (!user->noeqflag) {
201:     PetscCall(VecCreate(PETSC_COMM_WORLD, &user->ce)); /* a 1x1 vec for equality constraints */
202:     PetscCall(VecSetSizes(user->ce, neloc, user->ne));
203:     PetscCall(VecSetFromOptions(user->ce));
204:     PetscCall(VecSetUp(user->ce));
205:   }

207:   PetscCall(VecCreate(PETSC_COMM_WORLD, &user->ci)); /* a 2x1 vec for inequality constraints */
208:   PetscCall(VecSetSizes(user->ci, niloc, user->ni));
209:   PetscCall(VecSetFromOptions(user->ci));
210:   PetscCall(VecSetUp(user->ci));

212:   /* nexn & nixn matrices for equally and inequalty constraints */
213:   if (!user->noeqflag) {
214:     PetscCall(MatCreate(PETSC_COMM_WORLD, &user->Ae));
215:     PetscCall(MatSetSizes(user->Ae, neloc, nloc, user->ne, user->n));
216:     PetscCall(MatSetFromOptions(user->Ae));
217:     PetscCall(MatSetUp(user->Ae));
218:   }

220:   PetscCall(MatCreate(PETSC_COMM_WORLD, &user->Ai));
221:   PetscCall(MatSetSizes(user->Ai, niloc, nloc, user->ni, user->n));
222:   PetscCall(MatSetFromOptions(user->Ai));
223:   PetscCall(MatSetUp(user->Ai));

225:   PetscCall(MatCreate(PETSC_COMM_WORLD, &user->H));
226:   PetscCall(MatSetSizes(user->H, nloc, nloc, user->n, user->n));
227:   PetscCall(MatSetFromOptions(user->H));
228:   PetscCall(MatSetUp(user->H));
229:   PetscFunctionReturn(PETSC_SUCCESS);
230: }

232: PetscErrorCode DestroyProblem(AppCtx *user)
233: {
234:   PetscFunctionBegin;
235:   if (!user->noeqflag) PetscCall(MatDestroy(&user->Ae));
236:   PetscCall(MatDestroy(&user->Ai));
237:   PetscCall(MatDestroy(&user->H));

239:   PetscCall(VecDestroy(&user->x));
240:   if (!user->noeqflag) PetscCall(VecDestroy(&user->ce));
241:   PetscCall(VecDestroy(&user->ci));
242:   PetscCall(VecDestroy(&user->xl));
243:   PetscCall(VecDestroy(&user->xu));
244:   PetscCall(VecDestroy(&user->Xseq));
245:   PetscCall(VecScatterDestroy(&user->scat));
246:   PetscFunctionReturn(PETSC_SUCCESS);
247: }

249: /* Evaluate
250:    f(x) = (x0 - 2)^2 + (x1 - 2)^2 - 2*(x0 + x1)
251:    G = grad f = [2*(x0 - 2) - 2, 2*(x1 - 2) - 2] = 2*X - 6
252: */
253: PetscErrorCode FormFunctionGradient(Tao tao, Vec X, PetscReal *f, Vec G, PetscCtx ctx)
254: {
255:   const PetscScalar *x;
256:   MPI_Comm           comm;
257:   PetscMPIInt        rank;
258:   PetscReal          fin;
259:   AppCtx            *user = (AppCtx *)ctx;
260:   Vec                Xseq = user->Xseq;
261:   VecScatter         scat = user->scat;

263:   PetscFunctionBegin;
264:   /* f = (x0 - 2)^2 + (x1 - 2)^2 - 2*(x0 + x1) */
265:   PetscCall(PetscObjectGetComm((PetscObject)tao, &comm));
266:   PetscCallMPI(MPI_Comm_rank(comm, &rank));

268:   PetscCall(VecScatterBegin(scat, X, Xseq, INSERT_VALUES, SCATTER_FORWARD));
269:   PetscCall(VecScatterEnd(scat, X, Xseq, INSERT_VALUES, SCATTER_FORWARD));

271:   fin = 0.0;
272:   if (rank == 0) {
273:     PetscCall(VecGetArrayRead(Xseq, &x));
274:     fin = (x[0] - 2.0) * (x[0] - 2.0) + (x[1] - 2.0) * (x[1] - 2.0) - 2.0 * (x[0] + x[1]);
275:     PetscCall(VecRestoreArrayRead(Xseq, &x));
276:   }
277:   PetscCallMPI(MPIU_Allreduce(&fin, f, 1, MPIU_REAL, MPIU_SUM, comm));

279:   /* G = 2*X - 6 */
280:   PetscCall(VecSet(G, -6.0));
281:   PetscCall(VecAXPY(G, 2.0, X));
282:   PetscFunctionReturn(PETSC_SUCCESS);
283: }

285: /* Evaluate PDIPM Hessian, see Eqn(22) in http://doi.org/10.1049/gtd2.12708
286:    H = Wxx = fxx        + Jacobian(grad g^T*DE)   - Jacobian(grad h^T*DI)]
287:            = fxx        + Jacobin([2*x0; 1]de[0]) - Jacobian([2*x0, -2*x0; -1, 1][di[0] di[1]]^T)
288:            = [2 0; 0 2] + [2*de[0]  0;      0  0] - [2*di[0]-2*di[1], 0; 0, 0]
289:            = [ 2*(1+de[0]-di[0]+di[1]), 0;
290:                           0,            2]
291: */
292: PetscErrorCode FormPDIPMHessian(Tao tao, Vec x, Mat H, Mat Hpre, PetscCtx ctx)
293: {
294:   AppCtx            *user = (AppCtx *)ctx;
295:   Vec                DE, DI;
296:   const PetscScalar *de, *di;
297:   PetscInt           zero = 0, one = 1;
298:   PetscScalar        two = 2.0;
299:   PetscScalar        val = 0.0;
300:   Vec                Deseq, Diseq;
301:   VecScatter         Descat, Discat;
302:   PetscMPIInt        rank;
303:   MPI_Comm           comm;

305:   PetscFunctionBegin;
306:   PetscCall(TaoGetDualVariables(tao, &DE, &DI));

308:   PetscCall(PetscObjectGetComm((PetscObject)tao, &comm));
309:   PetscCallMPI(MPI_Comm_rank(comm, &rank));

311:   if (!user->noeqflag) {
312:     PetscCall(VecScatterCreateToZero(DE, &Descat, &Deseq));
313:     PetscCall(VecScatterBegin(Descat, DE, Deseq, INSERT_VALUES, SCATTER_FORWARD));
314:     PetscCall(VecScatterEnd(Descat, DE, Deseq, INSERT_VALUES, SCATTER_FORWARD));
315:   }
316:   PetscCall(VecScatterCreateToZero(DI, &Discat, &Diseq));
317:   PetscCall(VecScatterBegin(Discat, DI, Diseq, INSERT_VALUES, SCATTER_FORWARD));
318:   PetscCall(VecScatterEnd(Discat, DI, Diseq, INSERT_VALUES, SCATTER_FORWARD));

320:   if (rank == 0) {
321:     if (!user->noeqflag) PetscCall(VecGetArrayRead(Deseq, &de)); /* places equality constraint dual into array */
322:     PetscCall(VecGetArrayRead(Diseq, &di));                      /* places inequality constraint dual into array */

324:     if (!user->noeqflag) {
325:       val = 2.0 * (1 + de[0] - di[0] + di[1]);
326:       PetscCall(VecRestoreArrayRead(Deseq, &de));
327:       PetscCall(VecRestoreArrayRead(Diseq, &di));
328:     } else {
329:       val = 2.0 * (1 - di[0] + di[1]);
330:     }
331:     PetscCall(VecRestoreArrayRead(Diseq, &di));
332:     PetscCall(MatSetValues(H, 1, &zero, 1, &zero, &val, INSERT_VALUES));
333:     PetscCall(MatSetValues(H, 1, &one, 1, &one, &two, INSERT_VALUES));
334:   }
335:   PetscCall(MatAssemblyBegin(H, MAT_FINAL_ASSEMBLY));
336:   PetscCall(MatAssemblyEnd(H, MAT_FINAL_ASSEMBLY));
337:   if (!user->noeqflag) {
338:     PetscCall(VecScatterDestroy(&Descat));
339:     PetscCall(VecDestroy(&Deseq));
340:   }
341:   PetscCall(VecScatterDestroy(&Discat));
342:   PetscCall(VecDestroy(&Diseq));
343:   PetscFunctionReturn(PETSC_SUCCESS);
344: }

346: /* Evaluate
347:    h = [ x0^2 - x1;
348:          1 -(x0^2 - x1)]
349: */
350: PetscErrorCode FormInequalityConstraints(Tao tao, Vec X, Vec CI, PetscCtx ctx)
351: {
352:   const PetscScalar *x;
353:   PetscScalar        ci;
354:   MPI_Comm           comm;
355:   PetscMPIInt        rank;
356:   AppCtx            *user = (AppCtx *)ctx;
357:   Vec                Xseq = user->Xseq;
358:   VecScatter         scat = user->scat;

360:   PetscFunctionBegin;
361:   PetscCall(PetscObjectGetComm((PetscObject)tao, &comm));
362:   PetscCallMPI(MPI_Comm_rank(comm, &rank));

364:   PetscCall(VecScatterBegin(scat, X, Xseq, INSERT_VALUES, SCATTER_FORWARD));
365:   PetscCall(VecScatterEnd(scat, X, Xseq, INSERT_VALUES, SCATTER_FORWARD));

367:   if (rank == 0) {
368:     PetscCall(VecGetArrayRead(Xseq, &x));
369:     ci = x[0] * x[0] - x[1];
370:     PetscCall(VecSetValue(CI, 0, ci, INSERT_VALUES));
371:     ci = -x[0] * x[0] + x[1] + 1.0;
372:     PetscCall(VecSetValue(CI, 1, ci, INSERT_VALUES));
373:     PetscCall(VecRestoreArrayRead(Xseq, &x));
374:   }
375:   PetscCall(VecAssemblyBegin(CI));
376:   PetscCall(VecAssemblyEnd(CI));
377:   PetscFunctionReturn(PETSC_SUCCESS);
378: }

380: /* Evaluate
381:    g = [ x0^2 + x1 - 2]
382: */
383: PetscErrorCode FormEqualityConstraints(Tao tao, Vec X, Vec CE, PetscCtx ctx)
384: {
385:   const PetscScalar *x;
386:   PetscScalar        ce;
387:   MPI_Comm           comm;
388:   PetscMPIInt        rank;
389:   AppCtx            *user = (AppCtx *)ctx;
390:   Vec                Xseq = user->Xseq;
391:   VecScatter         scat = user->scat;

393:   PetscFunctionBegin;
394:   PetscCall(PetscObjectGetComm((PetscObject)tao, &comm));
395:   PetscCallMPI(MPI_Comm_rank(comm, &rank));

397:   PetscCall(VecScatterBegin(scat, X, Xseq, INSERT_VALUES, SCATTER_FORWARD));
398:   PetscCall(VecScatterEnd(scat, X, Xseq, INSERT_VALUES, SCATTER_FORWARD));

400:   if (rank == 0) {
401:     PetscCall(VecGetArrayRead(Xseq, &x));
402:     ce = x[0] * x[0] + x[1] - 2.0;
403:     PetscCall(VecSetValue(CE, 0, ce, INSERT_VALUES));
404:     PetscCall(VecRestoreArrayRead(Xseq, &x));
405:   }
406:   PetscCall(VecAssemblyBegin(CE));
407:   PetscCall(VecAssemblyEnd(CE));
408:   PetscFunctionReturn(PETSC_SUCCESS);
409: }

411: /*
412:   grad h = [  2*x0, -1;
413:              -2*x0,  1]
414: */
415: PetscErrorCode FormInequalityJacobian(Tao tao, Vec X, Mat JI, Mat JIpre, PetscCtx ctx)
416: {
417:   AppCtx            *user = (AppCtx *)ctx;
418:   PetscInt           zero = 0, one = 1, cols[2];
419:   PetscScalar        vals[2];
420:   const PetscScalar *x;
421:   Vec                Xseq = user->Xseq;
422:   VecScatter         scat = user->scat;
423:   MPI_Comm           comm;
424:   PetscMPIInt        rank;

426:   PetscFunctionBegin;
427:   PetscCall(PetscObjectGetComm((PetscObject)tao, &comm));
428:   PetscCallMPI(MPI_Comm_rank(comm, &rank));
429:   PetscCall(VecScatterBegin(scat, X, Xseq, INSERT_VALUES, SCATTER_FORWARD));
430:   PetscCall(VecScatterEnd(scat, X, Xseq, INSERT_VALUES, SCATTER_FORWARD));

432:   PetscCall(VecGetArrayRead(Xseq, &x));
433:   if (rank == 0) {
434:     cols[0] = 0;
435:     cols[1] = 1;
436:     vals[0] = 2 * x[0];
437:     vals[1] = -1.0;
438:     PetscCall(MatSetValues(JI, 1, &zero, 2, cols, vals, INSERT_VALUES));
439:     vals[0] = -2 * x[0];
440:     vals[1] = 1.0;
441:     PetscCall(MatSetValues(JI, 1, &one, 2, cols, vals, INSERT_VALUES));
442:   }
443:   PetscCall(VecRestoreArrayRead(Xseq, &x));
444:   PetscCall(MatAssemblyBegin(JI, MAT_FINAL_ASSEMBLY));
445:   PetscCall(MatAssemblyEnd(JI, MAT_FINAL_ASSEMBLY));
446:   PetscFunctionReturn(PETSC_SUCCESS);
447: }

449: /*
450:   grad g = [2*x0, 1.0]
451: */
452: PetscErrorCode FormEqualityJacobian(Tao tao, Vec X, Mat JE, Mat JEpre, PetscCtx ctx)
453: {
454:   PetscInt           zero = 0, cols[2];
455:   PetscScalar        vals[2];
456:   const PetscScalar *x;
457:   PetscMPIInt        rank;
458:   MPI_Comm           comm;

460:   PetscFunctionBegin;
461:   PetscCall(PetscObjectGetComm((PetscObject)tao, &comm));
462:   PetscCallMPI(MPI_Comm_rank(comm, &rank));

464:   if (rank == 0) {
465:     PetscCall(VecGetArrayRead(X, &x));
466:     cols[0] = 0;
467:     cols[1] = 1;
468:     vals[0] = 2 * x[0];
469:     vals[1] = 1.0;
470:     PetscCall(MatSetValues(JE, 1, &zero, 2, cols, vals, INSERT_VALUES));
471:     PetscCall(VecRestoreArrayRead(X, &x));
472:   }
473:   PetscCall(MatAssemblyBegin(JE, MAT_FINAL_ASSEMBLY));
474:   PetscCall(MatAssemblyEnd(JE, MAT_FINAL_ASSEMBLY));
475:   PetscFunctionReturn(PETSC_SUCCESS);
476: }

478: /*TEST

480:    build:
481:       requires: !complex !defined(PETSC_USE_CXX)

483:    test:
484:       args: -tao_converged_reason -tao_gatol 1.e-6 -tao_type pdipm -tao_pdipm_kkt_shift_pd
485:       requires: mumps !single
486:       filter: sed  -e "s/CONVERGED_GATOL iterations *[0-9]\{1,\}/CONVERGED_GATOL/g"

488:    test:
489:       suffix: pdipm_2
490:       requires: mumps
491:       nsize: 2
492:       args: -tao_converged_reason -tao_gatol 1.e-6 -tao_type pdipm
493:       filter: sed  -e "s/CONVERGED_GATOL iterations *[0-9]\{1,\}/CONVERGED_GATOL/g"

495:    test:
496:       suffix: 2
497:       args: -tao_converged_reason
498:       filter: sed  -e "s/CONVERGED_GATOL iterations *[0-9]\{1,\}/CONVERGED_GATOL/g"

500:    test:
501:       suffix: 3
502:       args: -tao_converged_reason -no_eq
503:       filter: sed  -e "s/CONVERGED_GATOL iterations *[0-9]\{1,\}/CONVERGED_GATOL/g"

505:    test:
506:       suffix: 4
507:       args: -tao_converged_reason -tao_almm_type classic
508:       requires: !single
509:       filter: sed  -e "s/CONVERGED_GATOL iterations *[0-9]\{1,\}/CONVERGED_GATOL/g"

511:    test:
512:       suffix: 5
513:       args: -tao_converged_reason -tao_almm_type classic -no_eq
514:       requires: !single
515:       filter: sed  -e "s/CONVERGED_GATOL iterations *[0-9]\{1,\}/CONVERGED_GATOL/g"

517:    test:
518:       suffix: 6
519:       args: -tao_converged_reason -tao_almm_subsolver_tao_type bqnktr
520:       filter: sed  -e "s/CONVERGED_GATOL iterations *[0-9]\{1,\}/CONVERGED_GATOL/g"

522:    test:
523:       suffix: 7
524:       args: -tao_converged_reason -tao_almm_subsolver_tao_type bncg
525:       filter: sed  -e "s/CONVERGED_GATOL iterations *[0-9]\{1,\}/CONVERGED_GATOL/g"

527:    test:
528:       suffix: 8
529:       nsize: 2
530:       args: -tao_converged_reason
531:       filter: sed  -e "s/CONVERGED_GATOL iterations *[0-9]\{1,\}/CONVERGED_GATOL/g"

533:    test:
534:       suffix: 9
535:       nsize: 2
536:       args: -tao_converged_reason -vec_type cuda -mat_type aijcusparse
537:       requires: cuda
538:       filter: sed  -e "s/CONVERGED_GATOL iterations *[0-9]\{1,\}/CONVERGED_GATOL/g"

540:    test:
541:       suffix: 10
542:       args: -tao_converged_reason -tao_almm_type classic -no_eq -no_bound
543:       requires: !single
544:       filter: sed  -e "s/CONVERGED_GATOL iterations *[0-9]\{1,\}/CONVERGED_GATOL/g"

546:    test:
547:       suffix: 11
548:       args: -tao_converged_reason -tao_almm_type classic -no_bound
549:       requires: !single
550:       filter: sed  -e "s/CONVERGED_GATOL iterations *[0-9]\{1,\}/CONVERGED_GATOL/g"

552:    test:
553:       suffix: 12
554:       args: -tao_converged_reason -tao_almm_type phr -no_eq -no_bound
555:       requires: !single
556:       filter: sed  -e "s/CONVERGED_GATOL iterations *[0-9]\{1,\}/CONVERGED_GATOL/g"

558:    test:
559:       suffix: 13
560:       args: -tao_converged_reason -tao_almm_type phr -no_bound
561:       requires: !single
562:       filter: sed  -e "s/CONVERGED_GATOL iterations *[0-9]\{1,\}/CONVERGED_GATOL/g"

564: TEST*/