xref: /petsc/src/ts/tutorials/power_grid/ex9adj.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Basic equation for generator stability analysis.\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown /*F
5c4762a1bSJed Brown 
6c4762a1bSJed Brown \begin{eqnarray}
7c4762a1bSJed Brown                  \frac{d \theta}{dt} = \omega_b (\omega - \omega_s)
8c4762a1bSJed Brown                  \frac{2 H}{\omega_s}\frac{d \omega}{dt} & = & P_m - P_max \sin(\theta) -D(\omega - \omega_s)\\
9c4762a1bSJed Brown \end{eqnarray}
10c4762a1bSJed Brown 
11c4762a1bSJed Brown   Ensemble of initial conditions
12c4762a1bSJed 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
13c4762a1bSJed Brown 
14c4762a1bSJed Brown   Fault at .1 seconds
15c4762a1bSJed 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
16c4762a1bSJed Brown 
17c4762a1bSJed Brown   Initial conditions same as when fault is ended
18c4762a1bSJed 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
19c4762a1bSJed Brown 
20c4762a1bSJed Brown F*/
21c4762a1bSJed Brown 
22c4762a1bSJed Brown /*
23c4762a1bSJed Brown    Include "petscts.h" so that we can use TS solvers.  Note that this
24c4762a1bSJed Brown    file automatically includes:
25c4762a1bSJed Brown      petscsys.h       - base PETSc routines   petscvec.h - vectors
26c4762a1bSJed Brown      petscmat.h - matrices
27c4762a1bSJed Brown      petscis.h     - index sets            petscksp.h - Krylov subspace methods
28c4762a1bSJed Brown      petscviewer.h - viewers               petscpc.h  - preconditioners
29c4762a1bSJed Brown      petscksp.h   - linear solvers
30c4762a1bSJed Brown */
31c4762a1bSJed Brown 
32c4762a1bSJed Brown #include <petscts.h>
33c4762a1bSJed Brown 
34c4762a1bSJed Brown typedef struct {
35c4762a1bSJed Brown   PetscScalar H, D, omega_b, omega_s, Pmax, Pm, E, V, X, u_s, c;
36c4762a1bSJed Brown   PetscInt    beta;
37c4762a1bSJed Brown   PetscReal   tf, tcl;
38c4762a1bSJed Brown } AppCtx;
39c4762a1bSJed Brown 
40*9371c9d4SSatish Balay PetscErrorCode PostStepFunction(TS ts) {
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));
51c4762a1bSJed Brown   PetscFunctionReturn(0);
52c4762a1bSJed Brown }
53c4762a1bSJed Brown 
54c4762a1bSJed Brown /*
55c4762a1bSJed Brown      Defines the ODE passed to the ODE solver
56c4762a1bSJed Brown */
57*9371c9d4SSatish Balay static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec U, Vec F, AppCtx *ctx) {
58c4762a1bSJed Brown   PetscScalar       *f, Pmax;
59c4762a1bSJed Brown   const PetscScalar *u;
60c4762a1bSJed Brown 
61c4762a1bSJed Brown   PetscFunctionBegin;
62c4762a1bSJed Brown   /*  The next three lines allow us to access the entries of the vectors directly */
639566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
649566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
65c4762a1bSJed 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 */
66c4762a1bSJed Brown   else Pmax = ctx->Pmax;
67c4762a1bSJed Brown 
68c4762a1bSJed Brown   f[0] = ctx->omega_b * (u[1] - ctx->omega_s);
69c4762a1bSJed Brown   f[1] = (-Pmax * PetscSinScalar(u[0]) - ctx->D * (u[1] - ctx->omega_s) + ctx->Pm) * ctx->omega_s / (2.0 * ctx->H);
70c4762a1bSJed Brown 
719566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
729566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
73c4762a1bSJed Brown   PetscFunctionReturn(0);
74c4762a1bSJed Brown }
75c4762a1bSJed Brown 
76c4762a1bSJed Brown /*
77c4762a1bSJed Brown      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
78c4762a1bSJed Brown */
79*9371c9d4SSatish Balay static PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec U, Mat A, Mat B, AppCtx *ctx) {
80c4762a1bSJed Brown   PetscInt           rowcol[] = {0, 1};
81c4762a1bSJed Brown   PetscScalar        J[2][2], Pmax;
82c4762a1bSJed Brown   const PetscScalar *u;
83c4762a1bSJed Brown 
84c4762a1bSJed Brown   PetscFunctionBegin;
859566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
86c4762a1bSJed 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 */
87c4762a1bSJed Brown   else Pmax = ctx->Pmax;
88c4762a1bSJed Brown 
89*9371c9d4SSatish Balay   J[0][0] = 0;
90*9371c9d4SSatish Balay   J[0][1] = ctx->omega_b;
91*9371c9d4SSatish Balay   J[1][1] = -ctx->D * ctx->omega_s / (2.0 * ctx->H);
92*9371c9d4SSatish Balay   J[1][0] = -Pmax * PetscCosScalar(u[0]) * ctx->omega_s / (2.0 * ctx->H);
93c4762a1bSJed Brown 
949566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES));
959566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
96c4762a1bSJed Brown 
979566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
989566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
99c4762a1bSJed Brown   if (A != B) {
1009566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
1019566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
102c4762a1bSJed Brown   }
103c4762a1bSJed Brown   PetscFunctionReturn(0);
104c4762a1bSJed Brown }
105c4762a1bSJed Brown 
106*9371c9d4SSatish Balay static PetscErrorCode RHSJacobianP(TS ts, PetscReal t, Vec X, Mat A, void *ctx0) {
107c4762a1bSJed Brown   PetscInt    row[] = {0, 1}, col[] = {0};
108c4762a1bSJed Brown   PetscScalar J[2][1];
109c4762a1bSJed Brown   AppCtx     *ctx = (AppCtx *)ctx0;
110c4762a1bSJed Brown 
111c4762a1bSJed Brown   PetscFunctionBeginUser;
112c4762a1bSJed Brown   J[0][0] = 0;
113c4762a1bSJed Brown   J[1][0] = ctx->omega_s / (2.0 * ctx->H);
1149566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 2, row, 1, col, &J[0][0], INSERT_VALUES));
1159566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1169566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
117c4762a1bSJed Brown   PetscFunctionReturn(0);
118c4762a1bSJed Brown }
119c4762a1bSJed Brown 
120*9371c9d4SSatish Balay static PetscErrorCode CostIntegrand(TS ts, PetscReal t, Vec U, Vec R, AppCtx *ctx) {
121c4762a1bSJed Brown   PetscScalar       *r;
122c4762a1bSJed Brown   const PetscScalar *u;
123c4762a1bSJed Brown 
124c4762a1bSJed Brown   PetscFunctionBegin;
1259566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1269566063dSJacob Faibussowitsch   PetscCall(VecGetArray(R, &r));
1272f613bf5SBarry Smith   r[0] = ctx->c * PetscPowScalarInt(PetscMax(0., u[0] - ctx->u_s), ctx->beta);
1289566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(R, &r));
1299566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
130c4762a1bSJed Brown   PetscFunctionReturn(0);
131c4762a1bSJed Brown }
132c4762a1bSJed Brown 
133*9371c9d4SSatish Balay static PetscErrorCode DRDUJacobianTranspose(TS ts, PetscReal t, Vec U, Mat DRDU, Mat B, AppCtx *ctx) {
134c4762a1bSJed Brown   PetscScalar        ru[1];
135c4762a1bSJed Brown   const PetscScalar *u;
136c4762a1bSJed Brown   PetscInt           row[] = {0}, col[] = {0};
137c4762a1bSJed Brown 
138c4762a1bSJed Brown   PetscFunctionBegin;
1399566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1402f613bf5SBarry Smith   ru[0] = ctx->c * ctx->beta * PetscPowScalarInt(PetscMax(0., u[0] - ctx->u_s), ctx->beta - 1);
1419566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1429566063dSJacob Faibussowitsch   PetscCall(MatSetValues(DRDU, 1, row, 1, col, ru, INSERT_VALUES));
1439566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(DRDU, MAT_FINAL_ASSEMBLY));
1449566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(DRDU, MAT_FINAL_ASSEMBLY));
145c4762a1bSJed Brown   PetscFunctionReturn(0);
146c4762a1bSJed Brown }
147c4762a1bSJed Brown 
148*9371c9d4SSatish Balay static PetscErrorCode DRDPJacobianTranspose(TS ts, PetscReal t, Vec U, Mat DRDP, AppCtx *ctx) {
149c4762a1bSJed Brown   PetscFunctionBegin;
1509566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(DRDP));
1519566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(DRDP, MAT_FINAL_ASSEMBLY));
1529566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(DRDP, MAT_FINAL_ASSEMBLY));
153c4762a1bSJed Brown   PetscFunctionReturn(0);
154c4762a1bSJed Brown }
155c4762a1bSJed Brown 
156*9371c9d4SSatish Balay PetscErrorCode ComputeSensiP(Vec lambda, Vec mu, AppCtx *ctx) {
157c4762a1bSJed Brown   PetscScalar        sensip;
158c4762a1bSJed Brown   const PetscScalar *x, *y;
159c4762a1bSJed Brown 
160c4762a1bSJed Brown   PetscFunctionBegin;
1619566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(lambda, &x));
1629566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(mu, &y));
163c4762a1bSJed Brown   sensip = 1. / PetscSqrtScalar(1. - (ctx->Pm / ctx->Pmax) * (ctx->Pm / ctx->Pmax)) / ctx->Pmax * x[0] + y[0];
1649566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n sensitivity wrt parameter pm: %.7f \n", (double)sensip));
1659566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(lambda, &x));
1669566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(mu, &y));
167c4762a1bSJed Brown   PetscFunctionReturn(0);
168c4762a1bSJed Brown }
169c4762a1bSJed Brown 
170*9371c9d4SSatish Balay int main(int argc, char **argv) {
171c4762a1bSJed Brown   TS           ts, quadts; /* ODE integrator */
172c4762a1bSJed Brown   Vec          U;          /* solution will be stored here */
173c4762a1bSJed Brown   Mat          A;          /* Jacobian matrix */
174c4762a1bSJed Brown   Mat          Jacp;       /* Jacobian matrix */
175c4762a1bSJed Brown   Mat          DRDU, DRDP;
176c4762a1bSJed Brown   PetscMPIInt  size;
177c4762a1bSJed Brown   PetscInt     n = 2;
178c4762a1bSJed Brown   AppCtx       ctx;
179c4762a1bSJed Brown   PetscScalar *u;
180c4762a1bSJed Brown   PetscReal    du[2]    = {0.0, 0.0};
181c4762a1bSJed Brown   PetscBool    ensemble = PETSC_FALSE, flg1, flg2;
182c4762a1bSJed Brown   PetscReal    ftime;
183c4762a1bSJed Brown   PetscInt     steps;
184c4762a1bSJed Brown   PetscScalar *x_ptr, *y_ptr;
185c4762a1bSJed Brown   Vec          lambda[1], q, mu[1];
186c4762a1bSJed Brown 
187c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
188c4762a1bSJed Brown      Initialize program
189c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
190327415f7SBarry Smith   PetscFunctionBeginUser;
1919566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
1929566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
1933c633725SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs");
194c4762a1bSJed Brown 
195c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
196c4762a1bSJed Brown     Create necessary matrix and vectors
197c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1989566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
1999566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE));
2009566063dSJacob Faibussowitsch   PetscCall(MatSetType(A, MATDENSE));
2019566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
2029566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
203c4762a1bSJed Brown 
2049566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &U, NULL));
205c4762a1bSJed Brown 
2069566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &Jacp));
2079566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(Jacp, PETSC_DECIDE, PETSC_DECIDE, 2, 1));
2089566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(Jacp));
2099566063dSJacob Faibussowitsch   PetscCall(MatSetUp(Jacp));
210c4762a1bSJed Brown 
2119566063dSJacob Faibussowitsch   PetscCall(MatCreateDense(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, 1, 1, NULL, &DRDP));
2129566063dSJacob Faibussowitsch   PetscCall(MatSetUp(DRDP));
2139566063dSJacob Faibussowitsch   PetscCall(MatCreateDense(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, 1, 2, NULL, &DRDU));
2149566063dSJacob Faibussowitsch   PetscCall(MatSetUp(DRDU));
215c4762a1bSJed Brown 
216c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
217c4762a1bSJed Brown     Set runtime options
218c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
219d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Swing equation options", "");
220c4762a1bSJed Brown   {
221c4762a1bSJed Brown     ctx.beta    = 2;
222c4762a1bSJed Brown     ctx.c       = 10000.0;
223c4762a1bSJed Brown     ctx.u_s     = 1.0;
224c4762a1bSJed Brown     ctx.omega_s = 1.0;
225c4762a1bSJed Brown     ctx.omega_b = 120.0 * PETSC_PI;
226c4762a1bSJed Brown     ctx.H       = 5.0;
2279566063dSJacob Faibussowitsch     PetscCall(PetscOptionsScalar("-Inertia", "", "", ctx.H, &ctx.H, NULL));
228c4762a1bSJed Brown     ctx.D = 5.0;
2299566063dSJacob Faibussowitsch     PetscCall(PetscOptionsScalar("-D", "", "", ctx.D, &ctx.D, NULL));
230c4762a1bSJed Brown     ctx.E    = 1.1378;
231c4762a1bSJed Brown     ctx.V    = 1.0;
232c4762a1bSJed Brown     ctx.X    = 0.545;
233c4762a1bSJed Brown     ctx.Pmax = ctx.E * ctx.V / ctx.X;
2349566063dSJacob Faibussowitsch     PetscCall(PetscOptionsScalar("-Pmax", "", "", ctx.Pmax, &ctx.Pmax, NULL));
235c4762a1bSJed Brown     ctx.Pm = 1.1;
2369566063dSJacob Faibussowitsch     PetscCall(PetscOptionsScalar("-Pm", "", "", ctx.Pm, &ctx.Pm, NULL));
237c4762a1bSJed Brown     ctx.tf  = 0.1;
238c4762a1bSJed Brown     ctx.tcl = 0.2;
2399566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-tf", "Time to start fault", "", ctx.tf, &ctx.tf, NULL));
2409566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-tcl", "Time to end fault", "", ctx.tcl, &ctx.tcl, NULL));
2419566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-ensemble", "Run ensemble of different initial conditions", "", ensemble, &ensemble, NULL));
242c4762a1bSJed Brown     if (ensemble) {
243c4762a1bSJed Brown       ctx.tf  = -1;
244c4762a1bSJed Brown       ctx.tcl = -1;
245c4762a1bSJed Brown     }
246c4762a1bSJed Brown 
2479566063dSJacob Faibussowitsch     PetscCall(VecGetArray(U, &u));
248c4762a1bSJed Brown     u[0] = PetscAsinScalar(ctx.Pm / ctx.Pmax);
249c4762a1bSJed Brown     u[1] = 1.0;
2509566063dSJacob Faibussowitsch     PetscCall(PetscOptionsRealArray("-u", "Initial solution", "", u, &n, &flg1));
251c4762a1bSJed Brown     n = 2;
2529566063dSJacob Faibussowitsch     PetscCall(PetscOptionsRealArray("-du", "Perturbation in initial solution", "", du, &n, &flg2));
253c4762a1bSJed Brown     u[0] += du[0];
254c4762a1bSJed Brown     u[1] += du[1];
2559566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(U, &u));
256c4762a1bSJed Brown     if (flg1 || flg2) {
257c4762a1bSJed Brown       ctx.tf  = -1;
258c4762a1bSJed Brown       ctx.tcl = -1;
259c4762a1bSJed Brown     }
260c4762a1bSJed Brown   }
261d0609cedSBarry Smith   PetscOptionsEnd();
262c4762a1bSJed Brown 
263c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
264c4762a1bSJed Brown      Create timestepping solver context
265c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2669566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
2679566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
2689566063dSJacob Faibussowitsch   PetscCall(TSSetEquationType(ts, TS_EQ_ODE_EXPLICIT)); /* less Jacobian evaluations when adjoint BEuler is used, otherwise no effect */
2699566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSRK));
2709566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ts, NULL, (TSRHSFunction)RHSFunction, &ctx));
2719566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobian(ts, A, A, (TSRHSJacobian)RHSJacobian, &ctx));
2729566063dSJacob Faibussowitsch   PetscCall(TSCreateQuadratureTS(ts, PETSC_TRUE, &quadts));
2739566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(quadts, NULL, (TSRHSFunction)CostIntegrand, &ctx));
2749566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobian(quadts, DRDU, DRDU, (TSRHSJacobian)DRDUJacobianTranspose, &ctx));
2759566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobianP(quadts, DRDP, (TSRHSJacobianP)DRDPJacobianTranspose, &ctx));
276c4762a1bSJed Brown 
277c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
278c4762a1bSJed Brown      Set initial conditions
279c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2809566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, U));
281c4762a1bSJed Brown 
282c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
283c4762a1bSJed Brown     Save trajectory of solution so that TSAdjointSolve() may be used
284c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2859566063dSJacob Faibussowitsch   PetscCall(TSSetSaveTrajectory(ts));
286c4762a1bSJed Brown 
2879566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &lambda[0], NULL));
288c4762a1bSJed Brown   /*   Set initial conditions for the adjoint integration */
2899566063dSJacob Faibussowitsch   PetscCall(VecGetArray(lambda[0], &y_ptr));
290*9371c9d4SSatish Balay   y_ptr[0] = 0.0;
291*9371c9d4SSatish Balay   y_ptr[1] = 0.0;
2929566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(lambda[0], &y_ptr));
293c4762a1bSJed Brown 
2949566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(Jacp, &mu[0], NULL));
2959566063dSJacob Faibussowitsch   PetscCall(VecGetArray(mu[0], &x_ptr));
296c4762a1bSJed Brown   x_ptr[0] = -1.0;
2979566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(mu[0], &x_ptr));
2989566063dSJacob Faibussowitsch   PetscCall(TSSetCostGradients(ts, 1, lambda, mu));
299c4762a1bSJed Brown 
300c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
301c4762a1bSJed Brown      Set solver options
302c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
3039566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, 10.0));
3049566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
3059566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, .01));
3069566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
307c4762a1bSJed Brown 
308c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
309c4762a1bSJed Brown      Solve nonlinear system
310c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
311c4762a1bSJed Brown   if (ensemble) {
312c4762a1bSJed Brown     for (du[1] = -2.5; du[1] <= .01; du[1] += .1) {
3139566063dSJacob Faibussowitsch       PetscCall(VecGetArray(U, &u));
314c4762a1bSJed Brown       u[0] = PetscAsinScalar(ctx.Pm / ctx.Pmax);
315c4762a1bSJed Brown       u[1] = ctx.omega_s;
316c4762a1bSJed Brown       u[0] += du[0];
317c4762a1bSJed Brown       u[1] += du[1];
3189566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(U, &u));
3199566063dSJacob Faibussowitsch       PetscCall(TSSetTimeStep(ts, .01));
3209566063dSJacob Faibussowitsch       PetscCall(TSSolve(ts, U));
321c4762a1bSJed Brown     }
322c4762a1bSJed Brown   } else {
3239566063dSJacob Faibussowitsch     PetscCall(TSSolve(ts, U));
324c4762a1bSJed Brown   }
3259566063dSJacob Faibussowitsch   PetscCall(VecView(U, PETSC_VIEWER_STDOUT_WORLD));
3269566063dSJacob Faibussowitsch   PetscCall(TSGetSolveTime(ts, &ftime));
3279566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts, &steps));
328c4762a1bSJed Brown 
329c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
330c4762a1bSJed Brown      Adjoint model starts here
331c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
332c4762a1bSJed Brown   /*   Set initial conditions for the adjoint integration */
3339566063dSJacob Faibussowitsch   PetscCall(VecGetArray(lambda[0], &y_ptr));
334*9371c9d4SSatish Balay   y_ptr[0] = 0.0;
335*9371c9d4SSatish Balay   y_ptr[1] = 0.0;
3369566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(lambda[0], &y_ptr));
337c4762a1bSJed Brown 
3389566063dSJacob Faibussowitsch   PetscCall(VecGetArray(mu[0], &x_ptr));
339c4762a1bSJed Brown   x_ptr[0] = -1.0;
3409566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(mu[0], &x_ptr));
341c4762a1bSJed Brown 
342c4762a1bSJed Brown   /*   Set RHS JacobianP */
3439566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobianP(ts, Jacp, RHSJacobianP, &ctx));
344c4762a1bSJed Brown 
3459566063dSJacob Faibussowitsch   PetscCall(TSAdjointSolve(ts));
346c4762a1bSJed Brown 
3479566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n sensitivity wrt initial conditions: d[Psi(tf)]/d[phi0]  d[Psi(tf)]/d[omega0]\n"));
3489566063dSJacob Faibussowitsch   PetscCall(VecView(lambda[0], PETSC_VIEWER_STDOUT_WORLD));
3499566063dSJacob Faibussowitsch   PetscCall(VecView(mu[0], PETSC_VIEWER_STDOUT_WORLD));
3509566063dSJacob Faibussowitsch   PetscCall(TSGetCostIntegral(ts, &q));
3519566063dSJacob Faibussowitsch   PetscCall(VecView(q, PETSC_VIEWER_STDOUT_WORLD));
3529566063dSJacob Faibussowitsch   PetscCall(VecGetArray(q, &x_ptr));
3539566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n cost function=%g\n", (double)(x_ptr[0] - ctx.Pm)));
3549566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(q, &x_ptr));
355c4762a1bSJed Brown 
3569566063dSJacob Faibussowitsch   PetscCall(ComputeSensiP(lambda[0], mu[0], &ctx));
357c4762a1bSJed Brown 
358c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
359c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
360c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
3619566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
3629566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Jacp));
3639566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&DRDU));
3649566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&DRDP));
3659566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&U));
3669566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&lambda[0]));
3679566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mu[0]));
3689566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
3699566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
370b122ec5aSJacob Faibussowitsch   return 0;
371c4762a1bSJed Brown }
372c4762a1bSJed Brown 
373c4762a1bSJed Brown /*TEST
374c4762a1bSJed Brown 
375c4762a1bSJed Brown    build:
376c4762a1bSJed Brown       requires: !complex
377c4762a1bSJed Brown 
378c4762a1bSJed Brown    test:
379c4762a1bSJed Brown       args: -viewer_binary_skip_info -ts_adapt_type none
380c4762a1bSJed Brown 
381c4762a1bSJed Brown TEST*/
382