xref: /petsc/src/ts/tutorials/power_grid/ex9opt.c (revision 1b37a2a7cc4a4fb30c3e967db1c694c0a1013f51)
1 static char help[] = "Basic equation for generator stability analysis.\n";
2 
3 /*F
4 
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}
9 
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
12 
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
15 
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
18 
19 F*/
20 
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 */
30 
31 #include <petsctao.h>
32 #include <petscts.h>
33 
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;
40 
41 PetscErrorCode FormFunction(Tao, Vec, PetscReal *, void *);
42 PetscErrorCode FormGradient(Tao, Vec, Vec, void *);
43 
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;
51 
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;
58 
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);
61 
62   PetscCall(VecRestoreArrayRead(U, &u));
63   PetscCall(VecRestoreArray(F, &f));
64   PetscFunctionReturn(PETSC_SUCCESS);
65 }
66 
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;
75 
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;
80 
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);
85 
86   PetscCall(MatSetValues(A, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES));
87   PetscCall(VecRestoreArrayRead(U, &u));
88 
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 }
97 
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;
103 
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 }
112 
113 static PetscErrorCode CostIntegrand(TS ts, PetscReal t, Vec U, Vec R, AppCtx *ctx)
114 {
115   PetscScalar       *r;
116   const PetscScalar *u;
117 
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 }
126 
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};
132 
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 }
142 
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 }
151 
152 PetscErrorCode ComputeSensiP(Vec lambda, Vec mu, AppCtx *ctx)
153 {
154   PetscScalar       *y, sensip;
155   const PetscScalar *x;
156 
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 }
166 
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;
183 
184   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
185      Initialize program
186      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
187   PetscFunctionBeginUser;
188   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
189   PetscFunctionBeginUser;
190   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
191   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
192 
193   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
194     Set runtime options
195     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
196   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Swing equation options", "");
197   {
198     ctx.beta    = 2;
199     ctx.c       = PetscRealConstant(10000.0);
200     ctx.u_s     = PetscRealConstant(1.0);
201     ctx.omega_s = PetscRealConstant(1.0);
202     ctx.omega_b = PetscRealConstant(120.0) * PETSC_PI;
203     ctx.H       = PetscRealConstant(5.0);
204     PetscCall(PetscOptionsScalar("-Inertia", "", "", ctx.H, &ctx.H, NULL));
205     ctx.D = PetscRealConstant(5.0);
206     PetscCall(PetscOptionsScalar("-D", "", "", ctx.D, &ctx.D, NULL));
207     ctx.E    = PetscRealConstant(1.1378);
208     ctx.V    = PetscRealConstant(1.0);
209     ctx.X    = PetscRealConstant(0.545);
210     ctx.Pmax = ctx.E * ctx.V / ctx.X;
211     PetscCall(PetscOptionsScalar("-Pmax", "", "", ctx.Pmax, &ctx.Pmax, NULL));
212     ctx.Pm = PetscRealConstant(1.0194);
213     PetscCall(PetscOptionsScalar("-Pm", "", "", ctx.Pm, &ctx.Pm, NULL));
214     ctx.tf  = PetscRealConstant(0.1);
215     ctx.tcl = PetscRealConstant(0.2);
216     PetscCall(PetscOptionsReal("-tf", "Time to start fault", "", ctx.tf, &ctx.tf, NULL));
217     PetscCall(PetscOptionsReal("-tcl", "Time to end fault", "", ctx.tcl, &ctx.tcl, NULL));
218   }
219   PetscOptionsEnd();
220 
221   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
222     Create necessary matrix and vectors
223     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
224   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
225   PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE));
226   PetscCall(MatSetType(A, MATDENSE));
227   PetscCall(MatSetFromOptions(A));
228   PetscCall(MatSetUp(A));
229 
230   PetscCall(MatCreateVecs(A, &U, NULL));
231 
232   PetscCall(MatCreate(PETSC_COMM_WORLD, &Jacp));
233   PetscCall(MatSetSizes(Jacp, PETSC_DECIDE, PETSC_DECIDE, 2, 1));
234   PetscCall(MatSetFromOptions(Jacp));
235   PetscCall(MatSetUp(Jacp));
236   PetscCall(MatCreateDense(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, 1, 1, NULL, &DRDP));
237   PetscCall(MatSetUp(DRDP));
238   PetscCall(MatCreateDense(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, 1, 2, NULL, &DRDU));
239   PetscCall(MatSetUp(DRDU));
240 
241   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
242      Create timestepping solver context
243      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
244   PetscCall(TSCreate(PETSC_COMM_WORLD, &ctx.ts));
245   PetscCall(TSSetProblemType(ctx.ts, TS_NONLINEAR));
246   PetscCall(TSSetEquationType(ctx.ts, TS_EQ_ODE_EXPLICIT)); /* less Jacobian evaluations when adjoint BEuler is used, otherwise no effect */
247   PetscCall(TSSetType(ctx.ts, TSRK));
248   PetscCall(TSSetRHSFunction(ctx.ts, NULL, (TSRHSFunctionFn *)RHSFunction, &ctx));
249   PetscCall(TSSetRHSJacobian(ctx.ts, A, A, (TSRHSJacobianFn *)RHSJacobian, &ctx));
250   PetscCall(TSSetExactFinalTime(ctx.ts, TS_EXACTFINALTIME_MATCHSTEP));
251 
252   PetscCall(MatCreateVecs(A, &lambda[0], NULL));
253   PetscCall(MatCreateVecs(Jacp, &mu[0], NULL));
254   PetscCall(TSSetCostGradients(ctx.ts, 1, lambda, mu));
255   PetscCall(TSSetRHSJacobianP(ctx.ts, Jacp, RHSJacobianP, &ctx));
256 
257   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
258      Set solver options
259    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
260   PetscCall(TSSetMaxTime(ctx.ts, PetscRealConstant(1.0)));
261   PetscCall(TSSetTimeStep(ctx.ts, PetscRealConstant(0.01)));
262   PetscCall(TSSetFromOptions(ctx.ts));
263 
264   PetscCall(TSGetTimeStep(ctx.ts, &ctx.dt)); /* save the stepsize */
265 
266   PetscCall(TSCreateQuadratureTS(ctx.ts, PETSC_TRUE, &quadts));
267   PetscCall(TSSetRHSFunction(quadts, NULL, (TSRHSFunctionFn *)CostIntegrand, &ctx));
268   PetscCall(TSSetRHSJacobian(quadts, DRDU, DRDU, (TSRHSJacobianFn *)DRDUJacobianTranspose, &ctx));
269   PetscCall(TSSetRHSJacobianP(quadts, DRDP, (TSRHSJacobianPFn *)DRDPJacobianTranspose, &ctx));
270   PetscCall(TSSetSolution(ctx.ts, U));
271 
272   /* Create TAO solver and set desired solution method */
273   PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao));
274   PetscCall(TaoSetType(tao, TAOBLMVM));
275 
276   /*
277      Optimization starts
278   */
279   /* Set initial solution guess */
280   PetscCall(VecCreateSeq(PETSC_COMM_WORLD, 1, &p));
281   PetscCall(VecGetArray(p, &x_ptr));
282   x_ptr[0] = ctx.Pm;
283   PetscCall(VecRestoreArray(p, &x_ptr));
284 
285   PetscCall(TaoSetSolution(tao, p));
286   /* Set routine for function and gradient evaluation */
287   PetscCall(TaoSetObjective(tao, FormFunction, (void *)&ctx));
288   PetscCall(TaoSetGradient(tao, NULL, FormGradient, (void *)&ctx));
289 
290   /* Set bounds for the optimization */
291   PetscCall(VecDuplicate(p, &lowerb));
292   PetscCall(VecDuplicate(p, &upperb));
293   PetscCall(VecGetArray(lowerb, &x_ptr));
294   x_ptr[0] = 0.;
295   PetscCall(VecRestoreArray(lowerb, &x_ptr));
296   PetscCall(VecGetArray(upperb, &x_ptr));
297   x_ptr[0] = PetscRealConstant(1.1);
298   PetscCall(VecRestoreArray(upperb, &x_ptr));
299   PetscCall(TaoSetVariableBounds(tao, lowerb, upperb));
300 
301   /* Check for any TAO command line options */
302   PetscCall(TaoSetFromOptions(tao));
303   PetscCall(TaoGetKSP(tao, &ksp));
304   if (ksp) {
305     PetscCall(KSPGetPC(ksp, &pc));
306     PetscCall(PCSetType(pc, PCNONE));
307   }
308 
309   /* SOLVE THE APPLICATION */
310   PetscCall(TaoSolve(tao));
311 
312   PetscCall(VecView(p, PETSC_VIEWER_STDOUT_WORLD));
313   /* Free TAO data structures */
314   PetscCall(TaoDestroy(&tao));
315   PetscCall(VecDestroy(&p));
316   PetscCall(VecDestroy(&lowerb));
317   PetscCall(VecDestroy(&upperb));
318 
319   PetscCall(TSDestroy(&ctx.ts));
320   PetscCall(VecDestroy(&U));
321   PetscCall(MatDestroy(&A));
322   PetscCall(MatDestroy(&Jacp));
323   PetscCall(MatDestroy(&DRDU));
324   PetscCall(MatDestroy(&DRDP));
325   PetscCall(VecDestroy(&lambda[0]));
326   PetscCall(VecDestroy(&mu[0]));
327   PetscCall(PetscFinalize());
328   return 0;
329 }
330 
331 /* ------------------------------------------------------------------ */
332 /*
333    FormFunction - Evaluates the function
334 
335    Input Parameters:
336    tao - the Tao context
337    X   - the input vector
338    ptr - optional user-defined context, as set by TaoSetObjectiveAndGradient()
339 
340    Output Parameters:
341    f   - the newly evaluated function
342 */
343 PetscErrorCode FormFunction(Tao tao, Vec P, PetscReal *f, void *ctx0)
344 {
345   AppCtx      *ctx = (AppCtx *)ctx0;
346   TS           ts  = ctx->ts;
347   Vec          U; /* solution will be stored here */
348   PetscScalar *u;
349   PetscScalar *x_ptr;
350   Vec          q;
351 
352   PetscFunctionBeginUser;
353   PetscCall(VecGetArrayRead(P, (const PetscScalar **)&x_ptr));
354   ctx->Pm = x_ptr[0];
355   PetscCall(VecRestoreArrayRead(P, (const PetscScalar **)&x_ptr));
356 
357   /* reset time */
358   PetscCall(TSSetTime(ts, 0.0));
359   /* reset step counter, this is critical for adjoint solver */
360   PetscCall(TSSetStepNumber(ts, 0));
361   /* reset step size, the step size becomes negative after TSAdjointSolve */
362   PetscCall(TSSetTimeStep(ts, ctx->dt));
363   /* reinitialize the integral value */
364   PetscCall(TSGetCostIntegral(ts, &q));
365   PetscCall(VecSet(q, 0.0));
366 
367   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
368      Set initial conditions
369    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
370   PetscCall(TSGetSolution(ts, &U));
371   PetscCall(VecGetArray(U, &u));
372   u[0] = PetscAsinScalar(ctx->Pm / ctx->Pmax);
373   u[1] = PetscRealConstant(1.0);
374   PetscCall(VecRestoreArray(U, &u));
375 
376   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
377      Solve nonlinear system
378      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
379   PetscCall(TSSolve(ts, U));
380   PetscCall(TSGetCostIntegral(ts, &q));
381   PetscCall(VecGetArray(q, &x_ptr));
382   *f = -ctx->Pm + x_ptr[0];
383   PetscCall(VecRestoreArray(q, &x_ptr));
384   PetscFunctionReturn(PETSC_SUCCESS);
385 }
386 
387 PetscErrorCode FormGradient(Tao tao, Vec P, Vec G, void *ctx0)
388 {
389   AppCtx      *ctx = (AppCtx *)ctx0;
390   TS           ts  = ctx->ts;
391   Vec          U; /* solution will be stored here */
392   PetscReal    ftime;
393   PetscInt     steps;
394   PetscScalar *u;
395   PetscScalar *x_ptr, *y_ptr;
396   Vec         *lambda, q, *mu;
397 
398   PetscFunctionBeginUser;
399   PetscCall(VecGetArrayRead(P, (const PetscScalar **)&x_ptr));
400   ctx->Pm = x_ptr[0];
401   PetscCall(VecRestoreArrayRead(P, (const PetscScalar **)&x_ptr));
402 
403   /* reset time */
404   PetscCall(TSSetTime(ts, 0.0));
405   /* reset step counter, this is critical for adjoint solver */
406   PetscCall(TSSetStepNumber(ts, 0));
407   /* reset step size, the step size becomes negative after TSAdjointSolve */
408   PetscCall(TSSetTimeStep(ts, ctx->dt));
409   /* reinitialize the integral value */
410   PetscCall(TSGetCostIntegral(ts, &q));
411   PetscCall(VecSet(q, 0.0));
412 
413   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
414      Set initial conditions
415    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
416   PetscCall(TSGetSolution(ts, &U));
417   PetscCall(VecGetArray(U, &u));
418   u[0] = PetscAsinScalar(ctx->Pm / ctx->Pmax);
419   u[1] = PetscRealConstant(1.0);
420   PetscCall(VecRestoreArray(U, &u));
421 
422   /* Set up to save trajectory before TSSetFromOptions() so that TSTrajectory options can be captured */
423   PetscCall(TSSetSaveTrajectory(ts));
424   PetscCall(TSSetFromOptions(ts));
425 
426   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
427      Solve nonlinear system
428      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
429   PetscCall(TSSolve(ts, U));
430 
431   PetscCall(TSGetSolveTime(ts, &ftime));
432   PetscCall(TSGetStepNumber(ts, &steps));
433 
434   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
435      Adjoint model starts here
436      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
437   PetscCall(TSGetCostGradients(ts, NULL, &lambda, &mu));
438   /*   Set initial conditions for the adjoint integration */
439   PetscCall(VecGetArray(lambda[0], &y_ptr));
440   y_ptr[0] = 0.0;
441   y_ptr[1] = 0.0;
442   PetscCall(VecRestoreArray(lambda[0], &y_ptr));
443   PetscCall(VecGetArray(mu[0], &x_ptr));
444   x_ptr[0] = PetscRealConstant(-1.0);
445   PetscCall(VecRestoreArray(mu[0], &x_ptr));
446 
447   PetscCall(TSAdjointSolve(ts));
448   PetscCall(TSGetCostIntegral(ts, &q));
449   PetscCall(ComputeSensiP(lambda[0], mu[0], ctx));
450   PetscCall(VecCopy(mu[0], G));
451   PetscFunctionReturn(PETSC_SUCCESS);
452 }
453 
454 /*TEST
455 
456    build:
457       requires: !complex !single
458 
459    test:
460       args: -viewer_binary_skip_info -ts_adapt_type none -tao_monitor -tao_gatol 0.0 -tao_grtol 1.e-3 -tao_converged_reason
461 
462    test:
463       suffix: 2
464       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 
466 TEST*/
467