Actual source code: ex9opt.c

  1: static char help[] = "Basic equation for generator stability analysis.\n";

  3: /*F

  5: \begin{eqnarray}
  6:                  \frac{d \theta}{dt} = \omega_b (\omega - \omega_s)
  7:                  \frac{2 H}{\omega_s}\frac{d \omega}{dt} & = & P_m - P_max \sin(\theta) -D(\omega - \omega_s)\\
  8: \end{eqnarray}

 10:   Ensemble of initial conditions
 11:    ./ex2 -ensemble -ts_monitor_draw_solution_phase -1,-3,3,3      -ts_adapt_dt_max .01  -ts_monitor -ts_type rosw -pc_type lu -ksp_type preonly

 13:   Fault at .1 seconds
 14:    ./ex2           -ts_monitor_draw_solution_phase .42,.95,.6,1.05 -ts_adapt_dt_max .01  -ts_monitor -ts_type rosw -pc_type lu -ksp_type preonly

 16:   Initial conditions same as when fault is ended
 17:    ./ex2 -u 0.496792,1.00932 -ts_monitor_draw_solution_phase .42,.95,.6,1.05  -ts_adapt_dt_max .01  -ts_monitor -ts_type rosw -pc_type lu -ksp_type preonly

 19: F*/

 21: /*
 22:    Include "petscts.h" so that we can use TS solvers.  Note that this
 23:    file automatically includes:
 24:      petscsys.h       - base PETSc routines   petscvec.h - vectors
 25:      petscmat.h - matrices
 26:      petscis.h     - index sets            petscksp.h - Krylov subspace methods
 27:      petscviewer.h - viewers               petscpc.h  - preconditioners
 28:      petscksp.h   - linear solvers
 29: */

 31: #include <petsctao.h>
 32: #include <petscts.h>

 34: typedef struct {
 35:   TS          ts;
 36:   PetscScalar H, D, omega_b, omega_s, Pmax, Pm, E, V, X, u_s, c;
 37:   PetscInt    beta;
 38:   PetscReal   tf, tcl, dt;
 39: } AppCtx;

 41: PetscErrorCode FormFunction(Tao, Vec, PetscReal *, void *);
 42: PetscErrorCode FormGradient(Tao, Vec, Vec, void *);

 44: /*
 45:      Defines the ODE passed to the ODE solver
 46: */
 47: static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec U, Vec F, AppCtx *ctx)
 48: {
 49:   PetscScalar       *f, Pmax;
 50:   const PetscScalar *u;

 52:   PetscFunctionBegin;
 53:   /*  The next three lines allow us to access the entries of the vectors directly */
 54:   PetscCall(VecGetArrayRead(U, &u));
 55:   PetscCall(VecGetArray(F, &f));
 56:   if ((t > ctx->tf) && (t < ctx->tcl)) Pmax = 0.0; /* A short-circuit on the generator terminal that drives the electrical power output (Pmax*sin(delta)) to 0 */
 57:   else Pmax = ctx->Pmax;

 59:   f[0] = ctx->omega_b * (u[1] - ctx->omega_s);
 60:   f[1] = (-Pmax * PetscSinScalar(u[0]) - ctx->D * (u[1] - ctx->omega_s) + ctx->Pm) * ctx->omega_s / (2.0 * ctx->H);

 62:   PetscCall(VecRestoreArrayRead(U, &u));
 63:   PetscCall(VecRestoreArray(F, &f));
 64:   PetscFunctionReturn(PETSC_SUCCESS);
 65: }

 67: /*
 68:      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
 69: */
 70: static PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec U, Mat A, Mat B, AppCtx *ctx)
 71: {
 72:   PetscInt           rowcol[] = {0, 1};
 73:   PetscScalar        J[2][2], Pmax;
 74:   const PetscScalar *u;

 76:   PetscFunctionBegin;
 77:   PetscCall(VecGetArrayRead(U, &u));
 78:   if ((t > ctx->tf) && (t < ctx->tcl)) Pmax = 0.0; /* A short-circuit on the generator terminal that drives the electrical power output (Pmax*sin(delta)) to 0 */
 79:   else Pmax = ctx->Pmax;

 81:   J[0][0] = 0;
 82:   J[0][1] = ctx->omega_b;
 83:   J[1][1] = -ctx->D * ctx->omega_s / (2.0 * ctx->H);
 84:   J[1][0] = -Pmax * PetscCosScalar(u[0]) * ctx->omega_s / (2.0 * ctx->H);

 86:   PetscCall(MatSetValues(A, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES));
 87:   PetscCall(VecRestoreArrayRead(U, &u));

 89:   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
 90:   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
 91:   if (A != B) {
 92:     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
 93:     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
 94:   }
 95:   PetscFunctionReturn(PETSC_SUCCESS);
 96: }

 98: static PetscErrorCode RHSJacobianP(TS ts, PetscReal t, Vec X, Mat A, void *ctx0)
 99: {
100:   PetscInt    row[] = {0, 1}, col[] = {0};
101:   PetscScalar J[2][1];
102:   AppCtx     *ctx = (AppCtx *)ctx0;

104:   PetscFunctionBeginUser;
105:   J[0][0] = 0;
106:   J[1][0] = ctx->omega_s / (2.0 * ctx->H);
107:   PetscCall(MatSetValues(A, 2, row, 1, col, &J[0][0], INSERT_VALUES));
108:   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
109:   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
110:   PetscFunctionReturn(PETSC_SUCCESS);
111: }

113: static PetscErrorCode CostIntegrand(TS ts, PetscReal t, Vec U, Vec R, AppCtx *ctx)
114: {
115:   PetscScalar       *r;
116:   const PetscScalar *u;

118:   PetscFunctionBegin;
119:   PetscCall(VecGetArrayRead(U, &u));
120:   PetscCall(VecGetArray(R, &r));
121:   r[0] = ctx->c * PetscPowScalarInt(PetscMax(0., u[0] - ctx->u_s), ctx->beta);
122:   PetscCall(VecRestoreArray(R, &r));
123:   PetscCall(VecRestoreArrayRead(U, &u));
124:   PetscFunctionReturn(PETSC_SUCCESS);
125: }

127: static PetscErrorCode DRDUJacobianTranspose(TS ts, PetscReal t, Vec U, Mat DRDU, Mat B, AppCtx *ctx)
128: {
129:   PetscScalar        ru[1];
130:   const PetscScalar *u;
131:   PetscInt           row[] = {0}, col[] = {0};

133:   PetscFunctionBegin;
134:   PetscCall(VecGetArrayRead(U, &u));
135:   ru[0] = ctx->c * ctx->beta * PetscPowScalarInt(PetscMax(0., u[0] - ctx->u_s), ctx->beta - 1);
136:   PetscCall(VecRestoreArrayRead(U, &u));
137:   PetscCall(MatSetValues(DRDU, 1, row, 1, col, ru, INSERT_VALUES));
138:   PetscCall(MatAssemblyBegin(DRDU, MAT_FINAL_ASSEMBLY));
139:   PetscCall(MatAssemblyEnd(DRDU, MAT_FINAL_ASSEMBLY));
140:   PetscFunctionReturn(PETSC_SUCCESS);
141: }

143: static PetscErrorCode DRDPJacobianTranspose(TS ts, PetscReal t, Vec U, Mat DRDP, AppCtx *ctx)
144: {
145:   PetscFunctionBegin;
146:   PetscCall(MatZeroEntries(DRDP));
147:   PetscCall(MatAssemblyBegin(DRDP, MAT_FINAL_ASSEMBLY));
148:   PetscCall(MatAssemblyEnd(DRDP, MAT_FINAL_ASSEMBLY));
149:   PetscFunctionReturn(PETSC_SUCCESS);
150: }

152: PetscErrorCode ComputeSensiP(Vec lambda, Vec mu, AppCtx *ctx)
153: {
154:   PetscScalar       *y, sensip;
155:   const PetscScalar *x;

157:   PetscFunctionBegin;
158:   PetscCall(VecGetArrayRead(lambda, &x));
159:   PetscCall(VecGetArray(mu, &y));
160:   sensip = 1. / PetscSqrtScalar(1. - (ctx->Pm / ctx->Pmax) * (ctx->Pm / ctx->Pmax)) / ctx->Pmax * x[0] + y[0];
161:   y[0]   = sensip;
162:   PetscCall(VecRestoreArray(mu, &y));
163:   PetscCall(VecRestoreArrayRead(lambda, &x));
164:   PetscFunctionReturn(PETSC_SUCCESS);
165: }

167: int main(int argc, char **argv)
168: {
169:   Vec          p;
170:   PetscScalar *x_ptr;
171:   PetscMPIInt  size;
172:   AppCtx       ctx;
173:   Vec          lowerb, upperb;
174:   Tao          tao;
175:   KSP          ksp;
176:   PC           pc;
177:   Vec          U, lambda[1], mu[1];
178:   Mat          A;    /* Jacobian matrix */
179:   Mat          Jacp; /* Jacobian matrix */
180:   Mat          DRDU, DRDP;
181:   PetscInt     n = 2;
182:   TS           quadts;

184:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
185:      Initialize program
186:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
187:   PetscFunctionBeginUser;
188:   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
189:   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
190:   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");

192:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
193:     Set runtime options
194:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
195:   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Swing equation options", "");
196:   {
197:     ctx.beta    = 2;
198:     ctx.c       = PetscRealConstant(10000.0);
199:     ctx.u_s     = PetscRealConstant(1.0);
200:     ctx.omega_s = PetscRealConstant(1.0);
201:     ctx.omega_b = PetscRealConstant(120.0) * PETSC_PI;
202:     ctx.H       = PetscRealConstant(5.0);
203:     PetscCall(PetscOptionsScalar("-Inertia", "", "", ctx.H, &ctx.H, NULL));
204:     ctx.D = PetscRealConstant(5.0);
205:     PetscCall(PetscOptionsScalar("-D", "", "", ctx.D, &ctx.D, NULL));
206:     ctx.E    = PetscRealConstant(1.1378);
207:     ctx.V    = PetscRealConstant(1.0);
208:     ctx.X    = PetscRealConstant(0.545);
209:     ctx.Pmax = ctx.E * ctx.V / ctx.X;
210:     PetscCall(PetscOptionsScalar("-Pmax", "", "", ctx.Pmax, &ctx.Pmax, NULL));
211:     ctx.Pm = PetscRealConstant(1.0194);
212:     PetscCall(PetscOptionsScalar("-Pm", "", "", ctx.Pm, &ctx.Pm, NULL));
213:     ctx.tf  = PetscRealConstant(0.1);
214:     ctx.tcl = PetscRealConstant(0.2);
215:     PetscCall(PetscOptionsReal("-tf", "Time to start fault", "", ctx.tf, &ctx.tf, NULL));
216:     PetscCall(PetscOptionsReal("-tcl", "Time to end fault", "", ctx.tcl, &ctx.tcl, NULL));
217:   }
218:   PetscOptionsEnd();

220:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
221:     Create necessary matrix and vectors
222:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
223:   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
224:   PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE));
225:   PetscCall(MatSetType(A, MATDENSE));
226:   PetscCall(MatSetFromOptions(A));
227:   PetscCall(MatSetUp(A));

229:   PetscCall(MatCreateVecs(A, &U, NULL));

231:   PetscCall(MatCreate(PETSC_COMM_WORLD, &Jacp));
232:   PetscCall(MatSetSizes(Jacp, PETSC_DECIDE, PETSC_DECIDE, 2, 1));
233:   PetscCall(MatSetFromOptions(Jacp));
234:   PetscCall(MatSetUp(Jacp));
235:   PetscCall(MatCreateDense(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, 1, 1, NULL, &DRDP));
236:   PetscCall(MatSetUp(DRDP));
237:   PetscCall(MatCreateDense(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, 1, 2, NULL, &DRDU));
238:   PetscCall(MatSetUp(DRDU));

240:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
241:      Create timestepping solver context
242:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
243:   PetscCall(TSCreate(PETSC_COMM_WORLD, &ctx.ts));
244:   PetscCall(TSSetProblemType(ctx.ts, TS_NONLINEAR));
245:   PetscCall(TSSetEquationType(ctx.ts, TS_EQ_ODE_EXPLICIT)); /* less Jacobian evaluations when adjoint BEuler is used, otherwise no effect */
246:   PetscCall(TSSetType(ctx.ts, TSRK));
247:   PetscCall(TSSetRHSFunction(ctx.ts, NULL, (TSRHSFunctionFn *)RHSFunction, &ctx));
248:   PetscCall(TSSetRHSJacobian(ctx.ts, A, A, (TSRHSJacobianFn *)RHSJacobian, &ctx));
249:   PetscCall(TSSetExactFinalTime(ctx.ts, TS_EXACTFINALTIME_MATCHSTEP));

251:   PetscCall(MatCreateVecs(A, &lambda[0], NULL));
252:   PetscCall(MatCreateVecs(Jacp, &mu[0], NULL));
253:   PetscCall(TSSetCostGradients(ctx.ts, 1, lambda, mu));
254:   PetscCall(TSSetRHSJacobianP(ctx.ts, Jacp, RHSJacobianP, &ctx));

256:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
257:      Set solver options
258:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
259:   PetscCall(TSSetMaxTime(ctx.ts, PetscRealConstant(1.0)));
260:   PetscCall(TSSetTimeStep(ctx.ts, PetscRealConstant(0.01)));
261:   PetscCall(TSSetFromOptions(ctx.ts));

263:   PetscCall(TSGetTimeStep(ctx.ts, &ctx.dt)); /* save the stepsize */

265:   PetscCall(TSCreateQuadratureTS(ctx.ts, PETSC_TRUE, &quadts));
266:   PetscCall(TSSetRHSFunction(quadts, NULL, (TSRHSFunctionFn *)CostIntegrand, &ctx));
267:   PetscCall(TSSetRHSJacobian(quadts, DRDU, DRDU, (TSRHSJacobianFn *)DRDUJacobianTranspose, &ctx));
268:   PetscCall(TSSetRHSJacobianP(quadts, DRDP, (TSRHSJacobianPFn *)DRDPJacobianTranspose, &ctx));
269:   PetscCall(TSSetSolution(ctx.ts, U));

271:   /* Create TAO solver and set desired solution method */
272:   PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao));
273:   PetscCall(TaoSetType(tao, TAOBLMVM));

275:   /*
276:      Optimization starts
277:   */
278:   /* Set initial solution guess */
279:   PetscCall(VecCreateSeq(PETSC_COMM_WORLD, 1, &p));
280:   PetscCall(VecGetArray(p, &x_ptr));
281:   x_ptr[0] = ctx.Pm;
282:   PetscCall(VecRestoreArray(p, &x_ptr));

284:   PetscCall(TaoSetSolution(tao, p));
285:   /* Set routine for function and gradient evaluation */
286:   PetscCall(TaoSetObjective(tao, FormFunction, (void *)&ctx));
287:   PetscCall(TaoSetGradient(tao, NULL, FormGradient, (void *)&ctx));

289:   /* Set bounds for the optimization */
290:   PetscCall(VecDuplicate(p, &lowerb));
291:   PetscCall(VecDuplicate(p, &upperb));
292:   PetscCall(VecGetArray(lowerb, &x_ptr));
293:   x_ptr[0] = 0.;
294:   PetscCall(VecRestoreArray(lowerb, &x_ptr));
295:   PetscCall(VecGetArray(upperb, &x_ptr));
296:   x_ptr[0] = PetscRealConstant(1.1);
297:   PetscCall(VecRestoreArray(upperb, &x_ptr));
298:   PetscCall(TaoSetVariableBounds(tao, lowerb, upperb));

300:   /* Check for any TAO command line options */
301:   PetscCall(TaoSetFromOptions(tao));
302:   PetscCall(TaoGetKSP(tao, &ksp));
303:   if (ksp) {
304:     PetscCall(KSPGetPC(ksp, &pc));
305:     PetscCall(PCSetType(pc, PCNONE));
306:   }

308:   /* SOLVE THE APPLICATION */
309:   PetscCall(TaoSolve(tao));

311:   PetscCall(VecView(p, PETSC_VIEWER_STDOUT_WORLD));
312:   /* Free TAO data structures */
313:   PetscCall(TaoDestroy(&tao));
314:   PetscCall(VecDestroy(&p));
315:   PetscCall(VecDestroy(&lowerb));
316:   PetscCall(VecDestroy(&upperb));

318:   PetscCall(TSDestroy(&ctx.ts));
319:   PetscCall(VecDestroy(&U));
320:   PetscCall(MatDestroy(&A));
321:   PetscCall(MatDestroy(&Jacp));
322:   PetscCall(MatDestroy(&DRDU));
323:   PetscCall(MatDestroy(&DRDP));
324:   PetscCall(VecDestroy(&lambda[0]));
325:   PetscCall(VecDestroy(&mu[0]));
326:   PetscCall(PetscFinalize());
327:   return 0;
328: }

330: /* ------------------------------------------------------------------ */
331: /*
332:    FormFunction - Evaluates the function

334:    Input Parameters:
335:    tao - the Tao context
336:    X   - the input vector
337:    ptr - optional user-defined context, as set by TaoSetObjectiveAndGradient()

339:    Output Parameters:
340:    f   - the newly evaluated function
341: */
342: PetscErrorCode FormFunction(Tao tao, Vec P, PetscReal *f, void *ctx0)
343: {
344:   AppCtx      *ctx = (AppCtx *)ctx0;
345:   TS           ts  = ctx->ts;
346:   Vec          U; /* solution will be stored here */
347:   PetscScalar *u;
348:   PetscScalar *x_ptr;
349:   Vec          q;

351:   PetscFunctionBeginUser;
352:   PetscCall(VecGetArrayRead(P, (const PetscScalar **)&x_ptr));
353:   ctx->Pm = x_ptr[0];
354:   PetscCall(VecRestoreArrayRead(P, (const PetscScalar **)&x_ptr));

356:   /* reset time */
357:   PetscCall(TSSetTime(ts, 0.0));
358:   /* reset step counter, this is critical for adjoint solver */
359:   PetscCall(TSSetStepNumber(ts, 0));
360:   /* reset step size, the step size becomes negative after TSAdjointSolve */
361:   PetscCall(TSSetTimeStep(ts, ctx->dt));
362:   /* reinitialize the integral value */
363:   PetscCall(TSGetCostIntegral(ts, &q));
364:   PetscCall(VecSet(q, 0.0));

366:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
367:      Set initial conditions
368:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
369:   PetscCall(TSGetSolution(ts, &U));
370:   PetscCall(VecGetArray(U, &u));
371:   u[0] = PetscAsinScalar(ctx->Pm / ctx->Pmax);
372:   u[1] = PetscRealConstant(1.0);
373:   PetscCall(VecRestoreArray(U, &u));

375:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
376:      Solve nonlinear system
377:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
378:   PetscCall(TSSolve(ts, U));
379:   PetscCall(TSGetCostIntegral(ts, &q));
380:   PetscCall(VecGetArray(q, &x_ptr));
381:   *f = -ctx->Pm + x_ptr[0];
382:   PetscCall(VecRestoreArray(q, &x_ptr));
383:   PetscFunctionReturn(PETSC_SUCCESS);
384: }

386: PetscErrorCode FormGradient(Tao tao, Vec P, Vec G, void *ctx0)
387: {
388:   AppCtx      *ctx = (AppCtx *)ctx0;
389:   TS           ts  = ctx->ts;
390:   Vec          U; /* solution will be stored here */
391:   PetscReal    ftime;
392:   PetscInt     steps;
393:   PetscScalar *u;
394:   PetscScalar *x_ptr, *y_ptr;
395:   Vec         *lambda, q, *mu;

397:   PetscFunctionBeginUser;
398:   PetscCall(VecGetArrayRead(P, (const PetscScalar **)&x_ptr));
399:   ctx->Pm = x_ptr[0];
400:   PetscCall(VecRestoreArrayRead(P, (const PetscScalar **)&x_ptr));

402:   /* reset time */
403:   PetscCall(TSSetTime(ts, 0.0));
404:   /* reset step counter, this is critical for adjoint solver */
405:   PetscCall(TSSetStepNumber(ts, 0));
406:   /* reset step size, the step size becomes negative after TSAdjointSolve */
407:   PetscCall(TSSetTimeStep(ts, ctx->dt));
408:   /* reinitialize the integral value */
409:   PetscCall(TSGetCostIntegral(ts, &q));
410:   PetscCall(VecSet(q, 0.0));

412:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
413:      Set initial conditions
414:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
415:   PetscCall(TSGetSolution(ts, &U));
416:   PetscCall(VecGetArray(U, &u));
417:   u[0] = PetscAsinScalar(ctx->Pm / ctx->Pmax);
418:   u[1] = PetscRealConstant(1.0);
419:   PetscCall(VecRestoreArray(U, &u));

421:   /* Set up to save trajectory before TSSetFromOptions() so that TSTrajectory options can be captured */
422:   PetscCall(TSSetSaveTrajectory(ts));
423:   PetscCall(TSSetFromOptions(ts));

425:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
426:      Solve nonlinear system
427:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
428:   PetscCall(TSSolve(ts, U));

430:   PetscCall(TSGetSolveTime(ts, &ftime));
431:   PetscCall(TSGetStepNumber(ts, &steps));

433:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
434:      Adjoint model starts here
435:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
436:   PetscCall(TSGetCostGradients(ts, NULL, &lambda, &mu));
437:   /*   Set initial conditions for the adjoint integration */
438:   PetscCall(VecGetArray(lambda[0], &y_ptr));
439:   y_ptr[0] = 0.0;
440:   y_ptr[1] = 0.0;
441:   PetscCall(VecRestoreArray(lambda[0], &y_ptr));
442:   PetscCall(VecGetArray(mu[0], &x_ptr));
443:   x_ptr[0] = PetscRealConstant(-1.0);
444:   PetscCall(VecRestoreArray(mu[0], &x_ptr));

446:   PetscCall(TSAdjointSolve(ts));
447:   PetscCall(TSGetCostIntegral(ts, &q));
448:   PetscCall(ComputeSensiP(lambda[0], mu[0], ctx));
449:   PetscCall(VecCopy(mu[0], G));
450:   PetscFunctionReturn(PETSC_SUCCESS);
451: }

453: /*TEST

455:    build:
456:       requires: !complex !single

458:    test:
459:       args: -viewer_binary_skip_info -ts_adapt_type none -tao_monitor -tao_gatol 0.0 -tao_grtol 1.e-3 -tao_converged_reason

461:    test:
462:       suffix: 2
463:       args: -viewer_binary_skip_info -ts_adapt_type none -tao_monitor -tao_gatol 0.0 -tao_grtol 1.e-3 -tao_converged_reason -tao_test_gradient

465: TEST*/