xref: /petsc/src/ts/tutorials/power_grid/ex9adj.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
1c4762a1bSJed Brown static char help[] = "Basic equation for generator stability analysis.\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown /*F
4c4762a1bSJed Brown 
5c4762a1bSJed Brown \begin{eqnarray}
6c4762a1bSJed Brown                  \frac{d \theta}{dt} = \omega_b (\omega - \omega_s)
7c4762a1bSJed Brown                  \frac{2 H}{\omega_s}\frac{d \omega}{dt} & = & P_m - P_max \sin(\theta) -D(\omega - \omega_s)\\
8c4762a1bSJed Brown \end{eqnarray}
9c4762a1bSJed Brown 
10c4762a1bSJed Brown   Ensemble of initial conditions
11c4762a1bSJed Brown    ./ex9 -ensemble -ts_monitor_draw_solution_phase -1,-3,3,3 -ts_adapt_dt_max .01 -ts_monitor -ts_type rk -pc_type lu -ksp_type preonly
12c4762a1bSJed Brown 
13c4762a1bSJed Brown   Fault at .1 seconds
14c4762a1bSJed Brown    ./ex9 -ts_monitor_draw_solution_phase .42,.95,.6,1.05 -ts_adapt_dt_max .01 -ts_monitor -ts_type rk -pc_type lu -ksp_type preonly
15c4762a1bSJed Brown 
16c4762a1bSJed Brown   Initial conditions same as when fault is ended
17c4762a1bSJed Brown    ./ex9 -u 0.496792,1.00932 -ts_monitor_draw_solution_phase .42,.95,.6,1.05 -ts_adapt_dt_max .01 -ts_monitor -ts_type rk -pc_type lu -ksp_type preonly
18c4762a1bSJed Brown 
19c4762a1bSJed Brown F*/
20c4762a1bSJed Brown 
21c4762a1bSJed Brown /*
22c4762a1bSJed Brown    Include "petscts.h" so that we can use TS solvers.  Note that this
23c4762a1bSJed Brown    file automatically includes:
24c4762a1bSJed Brown      petscsys.h       - base PETSc routines   petscvec.h - vectors
25c4762a1bSJed Brown      petscmat.h - matrices
26c4762a1bSJed Brown      petscis.h     - index sets            petscksp.h - Krylov subspace methods
27c4762a1bSJed Brown      petscviewer.h - viewers               petscpc.h  - preconditioners
28c4762a1bSJed Brown      petscksp.h   - linear solvers
29c4762a1bSJed Brown */
30c4762a1bSJed Brown 
31c4762a1bSJed Brown #include <petscts.h>
32c4762a1bSJed Brown 
33c4762a1bSJed Brown typedef struct {
34c4762a1bSJed Brown   PetscScalar H, D, omega_b, omega_s, Pmax, Pm, E, V, X, u_s, c;
35c4762a1bSJed Brown   PetscInt    beta;
36c4762a1bSJed Brown   PetscReal   tf, tcl;
37c4762a1bSJed Brown } AppCtx;
38c4762a1bSJed Brown 
PostStepFunction(TS ts)39d71ae5a4SJacob Faibussowitsch PetscErrorCode PostStepFunction(TS ts)
40d71ae5a4SJacob Faibussowitsch {
41c4762a1bSJed Brown   Vec                U;
42c4762a1bSJed Brown   PetscReal          t;
43c4762a1bSJed Brown   const PetscScalar *u;
44c4762a1bSJed Brown 
45c4762a1bSJed Brown   PetscFunctionBegin;
469566063dSJacob Faibussowitsch   PetscCall(TSGetTime(ts, &t));
479566063dSJacob Faibussowitsch   PetscCall(TSGetSolution(ts, &U));
489566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
499566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, "delta(%3.2f) = %8.7f\n", (double)t, (double)u[0]));
509566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
52c4762a1bSJed Brown }
53c4762a1bSJed Brown 
54c4762a1bSJed Brown /*
55c4762a1bSJed Brown      Defines the ODE passed to the ODE solver
56c4762a1bSJed Brown */
RHSFunction(TS ts,PetscReal t,Vec U,Vec F,AppCtx * ctx)57d71ae5a4SJacob Faibussowitsch static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec U, Vec F, AppCtx *ctx)
58d71ae5a4SJacob Faibussowitsch {
59c4762a1bSJed Brown   PetscScalar       *f, Pmax;
60c4762a1bSJed Brown   const PetscScalar *u;
61c4762a1bSJed Brown 
62c4762a1bSJed Brown   PetscFunctionBegin;
63c4762a1bSJed Brown   /*  The next three lines allow us to access the entries of the vectors directly */
649566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
659566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
66c4762a1bSJed Brown   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 */
67c4762a1bSJed Brown   else Pmax = ctx->Pmax;
68c4762a1bSJed Brown 
69c4762a1bSJed Brown   f[0] = ctx->omega_b * (u[1] - ctx->omega_s);
70c4762a1bSJed Brown   f[1] = (-Pmax * PetscSinScalar(u[0]) - ctx->D * (u[1] - ctx->omega_s) + ctx->Pm) * ctx->omega_s / (2.0 * ctx->H);
71c4762a1bSJed Brown 
729566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
739566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
75c4762a1bSJed Brown }
76c4762a1bSJed Brown 
77c4762a1bSJed Brown /*
78c4762a1bSJed Brown      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
79c4762a1bSJed Brown */
RHSJacobian(TS ts,PetscReal t,Vec U,Mat A,Mat B,AppCtx * ctx)80d71ae5a4SJacob Faibussowitsch static PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec U, Mat A, Mat B, AppCtx *ctx)
81d71ae5a4SJacob Faibussowitsch {
82c4762a1bSJed Brown   PetscInt           rowcol[] = {0, 1};
83c4762a1bSJed Brown   PetscScalar        J[2][2], Pmax;
84c4762a1bSJed Brown   const PetscScalar *u;
85c4762a1bSJed Brown 
86c4762a1bSJed Brown   PetscFunctionBegin;
879566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
88c4762a1bSJed Brown   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 */
89c4762a1bSJed Brown   else Pmax = ctx->Pmax;
90c4762a1bSJed Brown 
919371c9d4SSatish Balay   J[0][0] = 0;
929371c9d4SSatish Balay   J[0][1] = ctx->omega_b;
939371c9d4SSatish Balay   J[1][1] = -ctx->D * ctx->omega_s / (2.0 * ctx->H);
949371c9d4SSatish Balay   J[1][0] = -Pmax * PetscCosScalar(u[0]) * ctx->omega_s / (2.0 * ctx->H);
95c4762a1bSJed Brown 
969566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES));
979566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
98c4762a1bSJed Brown 
999566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1009566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
101c4762a1bSJed Brown   if (A != B) {
1029566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
1039566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
104c4762a1bSJed Brown   }
1053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
106c4762a1bSJed Brown }
107c4762a1bSJed Brown 
RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,PetscCtx ctx0)108*2a8381b2SBarry Smith static PetscErrorCode RHSJacobianP(TS ts, PetscReal t, Vec X, Mat A, PetscCtx ctx0)
109d71ae5a4SJacob Faibussowitsch {
110c4762a1bSJed Brown   PetscInt    row[] = {0, 1}, col[] = {0};
111c4762a1bSJed Brown   PetscScalar J[2][1];
112c4762a1bSJed Brown   AppCtx     *ctx = (AppCtx *)ctx0;
113c4762a1bSJed Brown 
114c4762a1bSJed Brown   PetscFunctionBeginUser;
115c4762a1bSJed Brown   J[0][0] = 0;
116c4762a1bSJed Brown   J[1][0] = ctx->omega_s / (2.0 * ctx->H);
1179566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 2, row, 1, col, &J[0][0], INSERT_VALUES));
1189566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1199566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
1203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
121c4762a1bSJed Brown }
122c4762a1bSJed Brown 
CostIntegrand(TS ts,PetscReal t,Vec U,Vec R,AppCtx * ctx)123d71ae5a4SJacob Faibussowitsch static PetscErrorCode CostIntegrand(TS ts, PetscReal t, Vec U, Vec R, AppCtx *ctx)
124d71ae5a4SJacob Faibussowitsch {
125c4762a1bSJed Brown   PetscScalar       *r;
126c4762a1bSJed Brown   const PetscScalar *u;
127c4762a1bSJed Brown 
128c4762a1bSJed Brown   PetscFunctionBegin;
1299566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1309566063dSJacob Faibussowitsch   PetscCall(VecGetArray(R, &r));
1312f613bf5SBarry Smith   r[0] = ctx->c * PetscPowScalarInt(PetscMax(0., u[0] - ctx->u_s), ctx->beta);
1329566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(R, &r));
1339566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
135c4762a1bSJed Brown }
136c4762a1bSJed Brown 
DRDUJacobianTranspose(TS ts,PetscReal t,Vec U,Mat DRDU,Mat B,AppCtx * ctx)137d71ae5a4SJacob Faibussowitsch static PetscErrorCode DRDUJacobianTranspose(TS ts, PetscReal t, Vec U, Mat DRDU, Mat B, AppCtx *ctx)
138d71ae5a4SJacob Faibussowitsch {
139c4762a1bSJed Brown   PetscScalar        ru[1];
140c4762a1bSJed Brown   const PetscScalar *u;
141c4762a1bSJed Brown   PetscInt           row[] = {0}, col[] = {0};
142c4762a1bSJed Brown 
143c4762a1bSJed Brown   PetscFunctionBegin;
1449566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1452f613bf5SBarry Smith   ru[0] = ctx->c * ctx->beta * PetscPowScalarInt(PetscMax(0., u[0] - ctx->u_s), ctx->beta - 1);
1469566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1479566063dSJacob Faibussowitsch   PetscCall(MatSetValues(DRDU, 1, row, 1, col, ru, INSERT_VALUES));
1489566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(DRDU, MAT_FINAL_ASSEMBLY));
1499566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(DRDU, MAT_FINAL_ASSEMBLY));
1503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
151c4762a1bSJed Brown }
152c4762a1bSJed Brown 
DRDPJacobianTranspose(TS ts,PetscReal t,Vec U,Mat DRDP,AppCtx * ctx)153d71ae5a4SJacob Faibussowitsch static PetscErrorCode DRDPJacobianTranspose(TS ts, PetscReal t, Vec U, Mat DRDP, AppCtx *ctx)
154d71ae5a4SJacob Faibussowitsch {
155c4762a1bSJed Brown   PetscFunctionBegin;
1569566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(DRDP));
1579566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(DRDP, MAT_FINAL_ASSEMBLY));
1589566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(DRDP, MAT_FINAL_ASSEMBLY));
1593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
160c4762a1bSJed Brown }
161c4762a1bSJed Brown 
ComputeSensiP(Vec lambda,Vec mu,AppCtx * ctx)162d71ae5a4SJacob Faibussowitsch PetscErrorCode ComputeSensiP(Vec lambda, Vec mu, AppCtx *ctx)
163d71ae5a4SJacob Faibussowitsch {
164c4762a1bSJed Brown   PetscScalar        sensip;
165c4762a1bSJed Brown   const PetscScalar *x, *y;
166c4762a1bSJed Brown 
167c4762a1bSJed Brown   PetscFunctionBegin;
1689566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(lambda, &x));
1699566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(mu, &y));
170c4762a1bSJed Brown   sensip = 1. / PetscSqrtScalar(1. - (ctx->Pm / ctx->Pmax) * (ctx->Pm / ctx->Pmax)) / ctx->Pmax * x[0] + y[0];
1719566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n sensitivity wrt parameter pm: %.7f \n", (double)sensip));
1729566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(lambda, &x));
1739566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(mu, &y));
1743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
175c4762a1bSJed Brown }
176c4762a1bSJed Brown 
main(int argc,char ** argv)177d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
178d71ae5a4SJacob Faibussowitsch {
179c4762a1bSJed Brown   TS           ts, quadts; /* ODE integrator */
180c4762a1bSJed Brown   Vec          U;          /* solution will be stored here */
181c4762a1bSJed Brown   Mat          A;          /* Jacobian matrix */
182c4762a1bSJed Brown   Mat          Jacp;       /* Jacobian matrix */
183c4762a1bSJed Brown   Mat          DRDU, DRDP;
184c4762a1bSJed Brown   PetscMPIInt  size;
185c4762a1bSJed Brown   PetscInt     n = 2;
186c4762a1bSJed Brown   AppCtx       ctx;
187c4762a1bSJed Brown   PetscScalar *u;
188c4762a1bSJed Brown   PetscReal    du[2]    = {0.0, 0.0};
189c4762a1bSJed Brown   PetscBool    ensemble = PETSC_FALSE, flg1, flg2;
190c4762a1bSJed Brown   PetscReal    ftime;
191c4762a1bSJed Brown   PetscInt     steps;
192c4762a1bSJed Brown   PetscScalar *x_ptr, *y_ptr;
193c4762a1bSJed Brown   Vec          lambda[1], q, mu[1];
194c4762a1bSJed Brown 
195c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
196c4762a1bSJed Brown      Initialize program
197c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
198327415f7SBarry Smith   PetscFunctionBeginUser;
199c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
2009566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
2013c633725SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs");
202c4762a1bSJed Brown 
203c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
204c4762a1bSJed Brown     Create necessary matrix and vectors
205c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2069566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
2079566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE));
2089566063dSJacob Faibussowitsch   PetscCall(MatSetType(A, MATDENSE));
2099566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
2109566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
211c4762a1bSJed Brown 
2129566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &U, NULL));
213c4762a1bSJed Brown 
2149566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &Jacp));
2159566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(Jacp, PETSC_DECIDE, PETSC_DECIDE, 2, 1));
2169566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(Jacp));
2179566063dSJacob Faibussowitsch   PetscCall(MatSetUp(Jacp));
218c4762a1bSJed Brown 
2199566063dSJacob Faibussowitsch   PetscCall(MatCreateDense(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, 1, 1, NULL, &DRDP));
2209566063dSJacob Faibussowitsch   PetscCall(MatSetUp(DRDP));
2219566063dSJacob Faibussowitsch   PetscCall(MatCreateDense(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, 1, 2, NULL, &DRDU));
2229566063dSJacob Faibussowitsch   PetscCall(MatSetUp(DRDU));
223c4762a1bSJed Brown 
224c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
225c4762a1bSJed Brown     Set runtime options
226c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
227d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Swing equation options", "");
228c4762a1bSJed Brown   {
229c4762a1bSJed Brown     ctx.beta    = 2;
230c4762a1bSJed Brown     ctx.c       = 10000.0;
231c4762a1bSJed Brown     ctx.u_s     = 1.0;
232c4762a1bSJed Brown     ctx.omega_s = 1.0;
233c4762a1bSJed Brown     ctx.omega_b = 120.0 * PETSC_PI;
234c4762a1bSJed Brown     ctx.H       = 5.0;
2359566063dSJacob Faibussowitsch     PetscCall(PetscOptionsScalar("-Inertia", "", "", ctx.H, &ctx.H, NULL));
236c4762a1bSJed Brown     ctx.D = 5.0;
2379566063dSJacob Faibussowitsch     PetscCall(PetscOptionsScalar("-D", "", "", ctx.D, &ctx.D, NULL));
238c4762a1bSJed Brown     ctx.E    = 1.1378;
239c4762a1bSJed Brown     ctx.V    = 1.0;
240c4762a1bSJed Brown     ctx.X    = 0.545;
241c4762a1bSJed Brown     ctx.Pmax = ctx.E * ctx.V / ctx.X;
2429566063dSJacob Faibussowitsch     PetscCall(PetscOptionsScalar("-Pmax", "", "", ctx.Pmax, &ctx.Pmax, NULL));
243c4762a1bSJed Brown     ctx.Pm = 1.1;
2449566063dSJacob Faibussowitsch     PetscCall(PetscOptionsScalar("-Pm", "", "", ctx.Pm, &ctx.Pm, NULL));
245c4762a1bSJed Brown     ctx.tf  = 0.1;
246c4762a1bSJed Brown     ctx.tcl = 0.2;
2479566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-tf", "Time to start fault", "", ctx.tf, &ctx.tf, NULL));
2489566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-tcl", "Time to end fault", "", ctx.tcl, &ctx.tcl, NULL));
2499566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-ensemble", "Run ensemble of different initial conditions", "", ensemble, &ensemble, NULL));
250c4762a1bSJed Brown     if (ensemble) {
251c4762a1bSJed Brown       ctx.tf  = -1;
252c4762a1bSJed Brown       ctx.tcl = -1;
253c4762a1bSJed Brown     }
254c4762a1bSJed Brown 
2559566063dSJacob Faibussowitsch     PetscCall(VecGetArray(U, &u));
256c4762a1bSJed Brown     u[0] = PetscAsinScalar(ctx.Pm / ctx.Pmax);
257c4762a1bSJed Brown     u[1] = 1.0;
2589566063dSJacob Faibussowitsch     PetscCall(PetscOptionsRealArray("-u", "Initial solution", "", u, &n, &flg1));
259c4762a1bSJed Brown     n = 2;
2609566063dSJacob Faibussowitsch     PetscCall(PetscOptionsRealArray("-du", "Perturbation in initial solution", "", du, &n, &flg2));
261c4762a1bSJed Brown     u[0] += du[0];
262c4762a1bSJed Brown     u[1] += du[1];
2639566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(U, &u));
264c4762a1bSJed Brown     if (flg1 || flg2) {
265c4762a1bSJed Brown       ctx.tf  = -1;
266c4762a1bSJed Brown       ctx.tcl = -1;
267c4762a1bSJed Brown     }
268c4762a1bSJed Brown   }
269d0609cedSBarry Smith   PetscOptionsEnd();
270c4762a1bSJed Brown 
271c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
272c4762a1bSJed Brown      Create timestepping solver context
273c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2749566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
2759566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
2769566063dSJacob Faibussowitsch   PetscCall(TSSetEquationType(ts, TS_EQ_ODE_EXPLICIT)); /* less Jacobian evaluations when adjoint BEuler is used, otherwise no effect */
2779566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSRK));
2788434afd1SBarry Smith   PetscCall(TSSetRHSFunction(ts, NULL, (TSRHSFunctionFn *)RHSFunction, &ctx));
2798434afd1SBarry Smith   PetscCall(TSSetRHSJacobian(ts, A, A, (TSRHSJacobianFn *)RHSJacobian, &ctx));
2809566063dSJacob Faibussowitsch   PetscCall(TSCreateQuadratureTS(ts, PETSC_TRUE, &quadts));
2818434afd1SBarry Smith   PetscCall(TSSetRHSFunction(quadts, NULL, (TSRHSFunctionFn *)CostIntegrand, &ctx));
2828434afd1SBarry Smith   PetscCall(TSSetRHSJacobian(quadts, DRDU, DRDU, (TSRHSJacobianFn *)DRDUJacobianTranspose, &ctx));
2838434afd1SBarry Smith   PetscCall(TSSetRHSJacobianP(quadts, DRDP, (TSRHSJacobianPFn *)DRDPJacobianTranspose, &ctx));
284fc905784SHong Zhang   PetscCall(TSSetCostGradients(ts, 1, lambda, mu));
285fc905784SHong Zhang   PetscCall(TSSetRHSJacobianP(ts, Jacp, RHSJacobianP, &ctx));
286c4762a1bSJed Brown 
287c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
288c4762a1bSJed Brown      Set initial conditions
289c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2909566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, U));
291c4762a1bSJed Brown 
292c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
293c4762a1bSJed Brown     Save trajectory of solution so that TSAdjointSolve() may be used
294c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2959566063dSJacob Faibussowitsch   PetscCall(TSSetSaveTrajectory(ts));
296c4762a1bSJed Brown 
2979566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &lambda[0], NULL));
298c4762a1bSJed Brown   /*   Set initial conditions for the adjoint integration */
2999566063dSJacob Faibussowitsch   PetscCall(VecGetArray(lambda[0], &y_ptr));
3009371c9d4SSatish Balay   y_ptr[0] = 0.0;
3019371c9d4SSatish Balay   y_ptr[1] = 0.0;
3029566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(lambda[0], &y_ptr));
303c4762a1bSJed Brown 
3049566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(Jacp, &mu[0], NULL));
3059566063dSJacob Faibussowitsch   PetscCall(VecGetArray(mu[0], &x_ptr));
306c4762a1bSJed Brown   x_ptr[0] = -1.0;
3079566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(mu[0], &x_ptr));
308c4762a1bSJed Brown 
309c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
310c4762a1bSJed Brown      Set solver options
311c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
3129566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, 10.0));
3139566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
3149566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, .01));
3159566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
316c4762a1bSJed Brown 
317c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
318c4762a1bSJed Brown      Solve nonlinear system
319c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
320c4762a1bSJed Brown   if (ensemble) {
321c4762a1bSJed Brown     for (du[1] = -2.5; du[1] <= .01; du[1] += .1) {
3229566063dSJacob Faibussowitsch       PetscCall(VecGetArray(U, &u));
323c4762a1bSJed Brown       u[0] = PetscAsinScalar(ctx.Pm / ctx.Pmax);
324c4762a1bSJed Brown       u[1] = ctx.omega_s;
325c4762a1bSJed Brown       u[0] += du[0];
326c4762a1bSJed Brown       u[1] += du[1];
3279566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(U, &u));
3289566063dSJacob Faibussowitsch       PetscCall(TSSetTimeStep(ts, .01));
3299566063dSJacob Faibussowitsch       PetscCall(TSSolve(ts, U));
330c4762a1bSJed Brown     }
331c4762a1bSJed Brown   } else {
3329566063dSJacob Faibussowitsch     PetscCall(TSSolve(ts, U));
333c4762a1bSJed Brown   }
3349566063dSJacob Faibussowitsch   PetscCall(VecView(U, PETSC_VIEWER_STDOUT_WORLD));
3359566063dSJacob Faibussowitsch   PetscCall(TSGetSolveTime(ts, &ftime));
3369566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts, &steps));
337c4762a1bSJed Brown 
338c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
339c4762a1bSJed Brown      Adjoint model starts here
340c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
341c4762a1bSJed Brown   /*   Set initial conditions for the adjoint integration */
3429566063dSJacob Faibussowitsch   PetscCall(VecGetArray(lambda[0], &y_ptr));
3439371c9d4SSatish Balay   y_ptr[0] = 0.0;
3449371c9d4SSatish Balay   y_ptr[1] = 0.0;
3459566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(lambda[0], &y_ptr));
346c4762a1bSJed Brown 
3479566063dSJacob Faibussowitsch   PetscCall(VecGetArray(mu[0], &x_ptr));
348c4762a1bSJed Brown   x_ptr[0] = -1.0;
3499566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(mu[0], &x_ptr));
350c4762a1bSJed Brown 
3519566063dSJacob Faibussowitsch   PetscCall(TSAdjointSolve(ts));
352c4762a1bSJed Brown 
3539566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n sensitivity wrt initial conditions: d[Psi(tf)]/d[phi0]  d[Psi(tf)]/d[omega0]\n"));
3549566063dSJacob Faibussowitsch   PetscCall(VecView(lambda[0], PETSC_VIEWER_STDOUT_WORLD));
3559566063dSJacob Faibussowitsch   PetscCall(VecView(mu[0], PETSC_VIEWER_STDOUT_WORLD));
3569566063dSJacob Faibussowitsch   PetscCall(TSGetCostIntegral(ts, &q));
3579566063dSJacob Faibussowitsch   PetscCall(VecView(q, PETSC_VIEWER_STDOUT_WORLD));
3589566063dSJacob Faibussowitsch   PetscCall(VecGetArray(q, &x_ptr));
3599566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n cost function=%g\n", (double)(x_ptr[0] - ctx.Pm)));
3609566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(q, &x_ptr));
361c4762a1bSJed Brown 
3629566063dSJacob Faibussowitsch   PetscCall(ComputeSensiP(lambda[0], mu[0], &ctx));
363c4762a1bSJed Brown 
364c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
365c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
366c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
3679566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
3689566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Jacp));
3699566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&DRDU));
3709566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&DRDP));
3719566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&U));
3729566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&lambda[0]));
3739566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mu[0]));
3749566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
3759566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
376b122ec5aSJacob Faibussowitsch   return 0;
377c4762a1bSJed Brown }
378c4762a1bSJed Brown 
379c4762a1bSJed Brown /*TEST
380c4762a1bSJed Brown 
381c4762a1bSJed Brown    build:
382c4762a1bSJed Brown       requires: !complex
383c4762a1bSJed Brown 
384c4762a1bSJed Brown    test:
385c4762a1bSJed Brown       args: -viewer_binary_skip_info -ts_adapt_type none
386c4762a1bSJed Brown 
387c4762a1bSJed Brown TEST*/
388