xref: /petsc/src/ts/tutorials/power_grid/ex9opt.c (revision 4e278199b78715991f5c71ebbd945c1489263e6c)
1 
2 static char help[] = "Basic equation for generator stability analysis.\n";
3 
4 /*F
5 
6 \begin{eqnarray}
7                  \frac{d \theta}{dt} = \omega_b (\omega - \omega_s)
8                  \frac{2 H}{\omega_s}\frac{d \omega}{dt} & = & P_m - P_max \sin(\theta) -D(\omega - \omega_s)\\
9 \end{eqnarray}
10 
11   Ensemble of initial conditions
12    ./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 
14   Fault at .1 seconds
15    ./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 
17   Initial conditions same as when fault is ended
18    ./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 
20 F*/
21 
22 /*
23    Include "petscts.h" so that we can use TS solvers.  Note that this
24    file automatically includes:
25      petscsys.h       - base PETSc routines   petscvec.h - vectors
26      petscmat.h - matrices
27      petscis.h     - index sets            petscksp.h - Krylov subspace methods
28      petscviewer.h - viewers               petscpc.h  - preconditioners
29      petscksp.h   - linear solvers
30 */
31 
32 #include <petsctao.h>
33 #include <petscts.h>
34 
35 typedef struct {
36   TS          ts;
37   PetscScalar H,D,omega_b,omega_s,Pmax,Pm,E,V,X,u_s,c;
38   PetscInt    beta;
39   PetscReal   tf,tcl,dt;
40 } AppCtx;
41 
42 PetscErrorCode FormFunction(Tao,Vec,PetscReal*,void*);
43 PetscErrorCode FormGradient(Tao,Vec,Vec,void*);
44 
45 /*
46      Defines the ODE passed to the ODE solver
47 */
48 static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec U,Vec F,AppCtx *ctx)
49 {
50   PetscErrorCode    ierr;
51   PetscScalar       *f,Pmax;
52   const PetscScalar *u;
53 
54   PetscFunctionBegin;
55   /*  The next three lines allow us to access the entries of the vectors directly */
56   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
57   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
58   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 */
59   else Pmax = ctx->Pmax;
60 
61   f[0] = ctx->omega_b*(u[1] - ctx->omega_s);
62   f[1] = (-Pmax*PetscSinScalar(u[0]) - ctx->D*(u[1] - ctx->omega_s) + ctx->Pm)*ctx->omega_s/(2.0*ctx->H);
63 
64   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
65   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
66   PetscFunctionReturn(0);
67 }
68 
69 /*
70      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
71 */
72 static PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec U,Mat A,Mat B,AppCtx *ctx)
73 {
74   PetscErrorCode    ierr;
75   PetscInt          rowcol[] = {0,1};
76   PetscScalar       J[2][2],Pmax;
77   const PetscScalar *u;
78 
79   PetscFunctionBegin;
80   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
81   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 */
82   else Pmax = ctx->Pmax;
83 
84   J[0][0] = 0;                                  J[0][1] = ctx->omega_b;
85   J[1][1] = -ctx->D*ctx->omega_s/(2.0*ctx->H);  J[1][0] = -Pmax*PetscCosScalar(u[0])*ctx->omega_s/(2.0*ctx->H);
86 
87   ierr    = MatSetValues(A,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);CHKERRQ(ierr);
88   ierr    = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
89 
90   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
91   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
92   if (A != B) {
93     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
94     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
95   }
96   PetscFunctionReturn(0);
97 }
98 
99 static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,void *ctx0)
100 {
101   PetscErrorCode ierr;
102   PetscInt       row[] = {0,1},col[]={0};
103   PetscScalar    J[2][1];
104   AppCtx         *ctx=(AppCtx*)ctx0;
105 
106   PetscFunctionBeginUser;
107   J[0][0] = 0;
108   J[1][0] = ctx->omega_s/(2.0*ctx->H);
109   ierr    = MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES);CHKERRQ(ierr);
110   ierr    = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
111   ierr    = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
112   PetscFunctionReturn(0);
113 }
114 
115 static PetscErrorCode CostIntegrand(TS ts,PetscReal t,Vec U,Vec R,AppCtx *ctx)
116 {
117   PetscErrorCode    ierr;
118   PetscScalar       *r;
119   const PetscScalar *u;
120 
121   PetscFunctionBegin;
122   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
123   ierr = VecGetArray(R,&r);CHKERRQ(ierr);
124   r[0] = ctx->c*PetscPowScalarInt(PetscMax(0., u[0]-ctx->u_s),ctx->beta);CHKERRQ(ierr);
125   ierr = VecRestoreArray(R,&r);CHKERRQ(ierr);
126   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
127   PetscFunctionReturn(0);
128 }
129 
130 static PetscErrorCode DRDUJacobianTranspose(TS ts,PetscReal t,Vec U,Mat DRDU,Mat B,AppCtx *ctx)
131 {
132   PetscErrorCode    ierr;
133   PetscScalar       ru[1];
134   const PetscScalar *u;
135   PetscInt          row[] = {0},col[] = {0};
136 
137   PetscFunctionBegin;
138   ierr  = VecGetArrayRead(U,&u);CHKERRQ(ierr);
139   ru[0] = ctx->c*ctx->beta*PetscPowScalarInt(PetscMax(0., u[0]-ctx->u_s),ctx->beta-1);CHKERRQ(ierr);
140   ierr  = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
141   ierr  = MatSetValues(DRDU,1,row,1,col,ru,INSERT_VALUES);CHKERRQ(ierr);
142   ierr  = MatAssemblyBegin(DRDU,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
143   ierr  = MatAssemblyEnd(DRDU,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
144   PetscFunctionReturn(0);
145 }
146 
147 static PetscErrorCode DRDPJacobianTranspose(TS ts,PetscReal t,Vec U,Mat DRDP,AppCtx *ctx)
148 {
149   PetscErrorCode ierr;
150 
151   PetscFunctionBegin;
152   ierr = MatZeroEntries(DRDP);CHKERRQ(ierr);
153   ierr = MatAssemblyBegin(DRDP,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
154   ierr = MatAssemblyEnd(DRDP,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
155   PetscFunctionReturn(0);
156 }
157 
158 PetscErrorCode ComputeSensiP(Vec lambda,Vec mu,AppCtx *ctx)
159 {
160   PetscErrorCode    ierr;
161   PetscScalar       *y,sensip;
162   const PetscScalar *x;
163 
164   PetscFunctionBegin;
165   ierr = VecGetArrayRead(lambda,&x);CHKERRQ(ierr);
166   ierr = VecGetArray(mu,&y);CHKERRQ(ierr);
167   sensip = 1./PetscSqrtScalar(1.-(ctx->Pm/ctx->Pmax)*(ctx->Pm/ctx->Pmax))/ctx->Pmax*x[0]+y[0];
168   y[0] = sensip;
169   ierr = VecRestoreArray(mu,&y);CHKERRQ(ierr);
170   ierr = VecRestoreArrayRead(lambda,&x);CHKERRQ(ierr);
171   PetscFunctionReturn(0);
172 }
173 
174 int main(int argc,char **argv)
175 {
176   Vec            p;
177   PetscScalar    *x_ptr;
178   PetscErrorCode ierr;
179   PetscMPIInt    size;
180   AppCtx         ctx;
181   Vec            lowerb,upperb;
182   Tao            tao;
183   KSP            ksp;
184   PC             pc;
185   Vec            U,lambda[1],mu[1];
186   Mat            A;             /* Jacobian matrix */
187   Mat            Jacp;          /* Jacobian matrix */
188   Mat            DRDU,DRDP;
189   PetscInt       n = 2;
190   TS             quadts;
191 
192   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
193      Initialize program
194      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
195   ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
196   PetscFunctionBeginUser;
197   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
198   if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
199 
200   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
201     Set runtime options
202     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
203   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Swing equation options","");CHKERRQ(ierr);
204   {
205     ctx.beta    = 2;
206     ctx.c       = PetscRealConstant(10000.0);
207     ctx.u_s     = PetscRealConstant(1.0);
208     ctx.omega_s = PetscRealConstant(1.0);
209     ctx.omega_b = PetscRealConstant(120.0)*PETSC_PI;
210     ctx.H       = PetscRealConstant(5.0);
211     ierr        = PetscOptionsScalar("-Inertia","","",ctx.H,&ctx.H,NULL);CHKERRQ(ierr);
212     ctx.D       = PetscRealConstant(5.0);
213     ierr        = PetscOptionsScalar("-D","","",ctx.D,&ctx.D,NULL);CHKERRQ(ierr);
214     ctx.E       = PetscRealConstant(1.1378);
215     ctx.V       = PetscRealConstant(1.0);
216     ctx.X       = PetscRealConstant(0.545);
217     ctx.Pmax    = ctx.E*ctx.V/ctx.X;
218     ierr        = PetscOptionsScalar("-Pmax","","",ctx.Pmax,&ctx.Pmax,NULL);CHKERRQ(ierr);
219     ctx.Pm      = PetscRealConstant(1.0194);
220     ierr        = PetscOptionsScalar("-Pm","","",ctx.Pm,&ctx.Pm,NULL);CHKERRQ(ierr);
221     ctx.tf      = PetscRealConstant(0.1);
222     ctx.tcl     = PetscRealConstant(0.2);
223     ierr        = PetscOptionsReal("-tf","Time to start fault","",ctx.tf,&ctx.tf,NULL);CHKERRQ(ierr);
224     ierr        = PetscOptionsReal("-tcl","Time to end fault","",ctx.tcl,&ctx.tcl,NULL);CHKERRQ(ierr);
225 
226   }
227   ierr = PetscOptionsEnd();CHKERRQ(ierr);
228 
229   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
230     Create necessary matrix and vectors
231     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
232   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
233   ierr = MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
234   ierr = MatSetType(A,MATDENSE);CHKERRQ(ierr);
235   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
236   ierr = MatSetUp(A);CHKERRQ(ierr);
237 
238   ierr = MatCreateVecs(A,&U,NULL);CHKERRQ(ierr);
239 
240   ierr = MatCreate(PETSC_COMM_WORLD,&Jacp);CHKERRQ(ierr);
241   ierr = MatSetSizes(Jacp,PETSC_DECIDE,PETSC_DECIDE,2,1);CHKERRQ(ierr);
242   ierr = MatSetFromOptions(Jacp);CHKERRQ(ierr);
243   ierr = MatSetUp(Jacp);CHKERRQ(ierr);
244   ierr = MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,&DRDP);CHKERRQ(ierr);
245   ierr = MatSetUp(DRDP);CHKERRQ(ierr);
246   ierr = MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,1,2,NULL,&DRDU);CHKERRQ(ierr);
247   ierr = MatSetUp(DRDU);CHKERRQ(ierr);
248 
249   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
250      Create timestepping solver context
251      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
252   ierr = TSCreate(PETSC_COMM_WORLD,&ctx.ts);CHKERRQ(ierr);
253   ierr = TSSetProblemType(ctx.ts,TS_NONLINEAR);CHKERRQ(ierr);
254   ierr = TSSetEquationType(ctx.ts,TS_EQ_ODE_EXPLICIT);CHKERRQ(ierr); /* less Jacobian evaluations when adjoint BEuler is used, otherwise no effect */
255   ierr = TSSetType(ctx.ts,TSRK);CHKERRQ(ierr);
256   ierr = TSSetRHSFunction(ctx.ts,NULL,(TSRHSFunction)RHSFunction,&ctx);CHKERRQ(ierr);
257   ierr = TSSetRHSJacobian(ctx.ts,A,A,(TSRHSJacobian)RHSJacobian,&ctx);CHKERRQ(ierr);
258   ierr = TSSetExactFinalTime(ctx.ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr);
259 
260   ierr = MatCreateVecs(A,&lambda[0],NULL);CHKERRQ(ierr);
261   ierr = MatCreateVecs(Jacp,&mu[0],NULL);CHKERRQ(ierr);
262   ierr = TSSetCostGradients(ctx.ts,1,lambda,mu);CHKERRQ(ierr);
263   ierr = TSSetRHSJacobianP(ctx.ts,Jacp,RHSJacobianP,&ctx);CHKERRQ(ierr);
264 
265   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
266      Set solver options
267    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
268   ierr = TSSetMaxTime(ctx.ts,PetscRealConstant(1.0));CHKERRQ(ierr);
269   ierr = TSSetTimeStep(ctx.ts,PetscRealConstant(0.01));CHKERRQ(ierr);
270   ierr = TSSetFromOptions(ctx.ts);CHKERRQ(ierr);
271 
272   ierr = TSGetTimeStep(ctx.ts,&ctx.dt);CHKERRQ(ierr); /* save the stepsize */
273 
274   ierr = TSCreateQuadratureTS(ctx.ts,PETSC_TRUE,&quadts);CHKERRQ(ierr);
275   ierr = TSSetRHSFunction(quadts,NULL,(TSRHSFunction)CostIntegrand,&ctx);CHKERRQ(ierr);
276   ierr = TSSetRHSJacobian(quadts,DRDU,DRDU,(TSRHSJacobian)DRDUJacobianTranspose,&ctx);CHKERRQ(ierr);
277   ierr = TSSetRHSJacobianP(quadts,DRDP,(TSRHSJacobianP)DRDPJacobianTranspose,&ctx);CHKERRQ(ierr);
278   ierr = TSSetSolution(ctx.ts,U);CHKERRQ(ierr);
279 
280   /* Create TAO solver and set desired solution method */
281   ierr = TaoCreate(PETSC_COMM_WORLD,&tao);CHKERRQ(ierr);
282   ierr = TaoSetType(tao,TAOBLMVM);CHKERRQ(ierr);
283 
284   /*
285      Optimization starts
286   */
287   /* Set initial solution guess */
288   ierr = VecCreateSeq(PETSC_COMM_WORLD,1,&p);CHKERRQ(ierr);
289   ierr = VecGetArray(p,&x_ptr);CHKERRQ(ierr);
290   x_ptr[0]   = ctx.Pm;
291   ierr = VecRestoreArray(p,&x_ptr);CHKERRQ(ierr);
292 
293   ierr = TaoSetInitialVector(tao,p);CHKERRQ(ierr);
294   /* Set routine for function and gradient evaluation */
295   ierr = TaoSetObjectiveRoutine(tao,FormFunction,(void *)&ctx);CHKERRQ(ierr);
296   ierr = TaoSetGradientRoutine(tao,FormGradient,(void *)&ctx);CHKERRQ(ierr);
297 
298   /* Set bounds for the optimization */
299   ierr = VecDuplicate(p,&lowerb);CHKERRQ(ierr);
300   ierr = VecDuplicate(p,&upperb);CHKERRQ(ierr);
301   ierr = VecGetArray(lowerb,&x_ptr);CHKERRQ(ierr);
302   x_ptr[0] = 0.;
303   ierr = VecRestoreArray(lowerb,&x_ptr);CHKERRQ(ierr);
304   ierr = VecGetArray(upperb,&x_ptr);CHKERRQ(ierr);
305   x_ptr[0] = PetscRealConstant(1.1);
306   ierr = VecRestoreArray(upperb,&x_ptr);CHKERRQ(ierr);
307   ierr = TaoSetVariableBounds(tao,lowerb,upperb);CHKERRQ(ierr);
308 
309   /* Check for any TAO command line options */
310   ierr = TaoSetFromOptions(tao);CHKERRQ(ierr);
311   ierr = TaoGetKSP(tao,&ksp);CHKERRQ(ierr);
312   if (ksp) {
313     ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
314     ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
315   }
316 
317   /* SOLVE THE APPLICATION */
318   ierr = TaoSolve(tao);CHKERRQ(ierr);
319 
320   ierr = VecView(p,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
321   /* Free TAO data structures */
322   ierr = TaoDestroy(&tao);CHKERRQ(ierr);
323   ierr = VecDestroy(&p);CHKERRQ(ierr);
324   ierr = VecDestroy(&lowerb);CHKERRQ(ierr);
325   ierr = VecDestroy(&upperb);CHKERRQ(ierr);
326 
327   ierr = TSDestroy(&ctx.ts);CHKERRQ(ierr);
328   ierr = VecDestroy(&U);CHKERRQ(ierr);
329   ierr = MatDestroy(&A);CHKERRQ(ierr);
330   ierr = MatDestroy(&Jacp);CHKERRQ(ierr);
331   ierr = MatDestroy(&DRDU);CHKERRQ(ierr);
332   ierr = MatDestroy(&DRDP);CHKERRQ(ierr);
333   ierr = VecDestroy(&lambda[0]);CHKERRQ(ierr);
334   ierr = VecDestroy(&mu[0]);CHKERRQ(ierr);
335   ierr = PetscFinalize();
336   return ierr;
337 }
338 
339 /* ------------------------------------------------------------------ */
340 /*
341    FormFunction - Evaluates the function
342 
343    Input Parameters:
344    tao - the Tao context
345    X   - the input vector
346    ptr - optional user-defined context, as set by TaoSetObjectiveAndGradientRoutine()
347 
348    Output Parameters:
349    f   - the newly evaluated function
350 */
351 PetscErrorCode FormFunction(Tao tao,Vec P,PetscReal *f,void *ctx0)
352 {
353   AppCtx         *ctx = (AppCtx*)ctx0;
354   TS             ts = ctx->ts;
355   Vec            U;             /* solution will be stored here */
356   PetscErrorCode ierr;
357   PetscScalar    *u;
358   PetscScalar    *x_ptr;
359   Vec            q;
360 
361   ierr = VecGetArrayRead(P,(const PetscScalar**)&x_ptr);CHKERRQ(ierr);
362   ctx->Pm = x_ptr[0];
363   ierr = VecRestoreArrayRead(P,(const PetscScalar**)&x_ptr);CHKERRQ(ierr);
364 
365   /* reset time */
366   ierr = TSSetTime(ts,0.0);CHKERRQ(ierr);
367   /* reset step counter, this is critical for adjoint solver */
368   ierr = TSSetStepNumber(ts,0);CHKERRQ(ierr);
369   /* reset step size, the step size becomes negative after TSAdjointSolve */
370   ierr = TSSetTimeStep(ts,ctx->dt);CHKERRQ(ierr);
371   /* reinitialize the integral value */
372   ierr = TSGetCostIntegral(ts,&q);CHKERRQ(ierr);
373   ierr = VecSet(q,0.0);CHKERRQ(ierr);
374 
375   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
376      Set initial conditions
377    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
378   ierr = TSGetSolution(ts,&U);CHKERRQ(ierr);
379   ierr = VecGetArray(U,&u);CHKERRQ(ierr);
380   u[0] = PetscAsinScalar(ctx->Pm/ctx->Pmax);
381   u[1] = PetscRealConstant(1.0);
382   ierr = VecRestoreArray(U,&u);CHKERRQ(ierr);
383 
384   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
385      Solve nonlinear system
386      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
387   ierr = TSSolve(ts,U);CHKERRQ(ierr);
388   ierr = TSGetCostIntegral(ts,&q);CHKERRQ(ierr);
389   ierr = VecGetArray(q,&x_ptr);CHKERRQ(ierr);
390   *f   = -ctx->Pm + x_ptr[0];
391   ierr = VecRestoreArray(q,&x_ptr);CHKERRQ(ierr);
392   return 0;
393 }
394 
395 PetscErrorCode FormGradient(Tao tao,Vec P,Vec G,void *ctx0)
396 {
397   AppCtx         *ctx = (AppCtx*)ctx0;
398   TS             ts = ctx->ts;
399   Vec            U;             /* solution will be stored here */
400   PetscErrorCode ierr;
401   PetscReal      ftime;
402   PetscInt       steps;
403   PetscScalar    *u;
404   PetscScalar    *x_ptr,*y_ptr;
405   Vec            *lambda,q,*mu;
406 
407   ierr = VecGetArrayRead(P,(const PetscScalar**)&x_ptr);CHKERRQ(ierr);
408   ctx->Pm = x_ptr[0];
409   ierr = VecRestoreArrayRead(P,(const PetscScalar**)&x_ptr);CHKERRQ(ierr);
410 
411   /* reset time */
412   ierr = TSSetTime(ts,0.0);CHKERRQ(ierr);
413   /* reset step counter, this is critical for adjoint solver */
414   ierr = TSSetStepNumber(ts,0);CHKERRQ(ierr);
415   /* reset step size, the step size becomes negative after TSAdjointSolve */
416   ierr = TSSetTimeStep(ts,ctx->dt);CHKERRQ(ierr);
417   /* reinitialize the integral value */
418   ierr = TSGetCostIntegral(ts,&q);CHKERRQ(ierr);
419   ierr = VecSet(q,0.0);CHKERRQ(ierr);
420 
421   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
422      Set initial conditions
423    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
424   ierr = TSGetSolution(ts,&U);CHKERRQ(ierr);
425   ierr = VecGetArray(U,&u);CHKERRQ(ierr);
426   u[0] = PetscAsinScalar(ctx->Pm/ctx->Pmax);
427   u[1] = PetscRealConstant(1.0);
428   ierr = VecRestoreArray(U,&u);CHKERRQ(ierr);
429 
430   /* Set up to save trajectory before TSSetFromOptions() so that TSTrajectory options can be captured */
431   ierr = TSSetSaveTrajectory(ts);CHKERRQ(ierr);
432   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
433 
434   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
435      Solve nonlinear system
436      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
437   ierr = TSSolve(ts,U);CHKERRQ(ierr);
438 
439   ierr = TSGetSolveTime(ts,&ftime);CHKERRQ(ierr);
440   ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr);
441 
442   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
443      Adjoint model starts here
444      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
445   ierr = TSGetCostGradients(ts,NULL,&lambda,&mu);CHKERRQ(ierr);
446   /*   Set initial conditions for the adjoint integration */
447   ierr = VecGetArray(lambda[0],&y_ptr);CHKERRQ(ierr);
448   y_ptr[0] = 0.0; y_ptr[1] = 0.0;
449   ierr = VecRestoreArray(lambda[0],&y_ptr);CHKERRQ(ierr);
450   ierr = VecGetArray(mu[0],&x_ptr);CHKERRQ(ierr);
451   x_ptr[0] = PetscRealConstant(-1.0);
452   ierr = VecRestoreArray(mu[0],&x_ptr);CHKERRQ(ierr);
453 
454   ierr = TSAdjointSolve(ts);CHKERRQ(ierr);
455   ierr = TSGetCostIntegral(ts,&q);CHKERRQ(ierr);
456   ierr = ComputeSensiP(lambda[0],mu[0],ctx);CHKERRQ(ierr);
457   ierr = VecCopy(mu[0],G);CHKERRQ(ierr);
458   return 0;
459 }
460 
461 /*TEST
462 
463    build:
464       requires: !complex
465 
466    test:
467       args: -viewer_binary_skip_info -ts_adapt_type none -tao_monitor -tao_gatol 0.0 -tao_grtol 1.e-3 -tao_converged_reason
468 
469    test:
470       suffix: 2
471       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
472 
473 TEST*/
474