1c4762a1bSJed Brown static char help[] = "Adjoint and tangent linear sensitivity analysis of the 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 F*/
11c4762a1bSJed Brown
12c4762a1bSJed Brown /*
13c4762a1bSJed Brown This code demonstrate the sensitivity analysis interface to a system of ordinary differential equations with discontinuities.
14c4762a1bSJed Brown It computes the sensitivities of an integral cost function
15c4762a1bSJed Brown \int c*max(0,\theta(t)-u_s)^beta dt
16c4762a1bSJed Brown w.r.t. initial conditions and the parameter P_m.
17c4762a1bSJed Brown Backward Euler method is used for time integration.
18c4762a1bSJed Brown The discontinuities are detected with TSEvent.
19c4762a1bSJed Brown */
20c4762a1bSJed Brown
21c4762a1bSJed Brown #include <petscts.h>
22c4762a1bSJed Brown #include "ex3.h"
23c4762a1bSJed Brown
main(int argc,char ** argv)24d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
25d71ae5a4SJacob Faibussowitsch {
26c4762a1bSJed Brown TS ts, quadts; /* ODE integrator */
27c4762a1bSJed Brown Vec U; /* solution will be stored here */
28c4762a1bSJed Brown PetscMPIInt size;
29c4762a1bSJed Brown PetscInt n = 2;
30c4762a1bSJed Brown AppCtx ctx;
31c4762a1bSJed Brown PetscScalar *u;
32c4762a1bSJed Brown PetscReal du[2] = {0.0, 0.0};
33c4762a1bSJed Brown PetscBool ensemble = PETSC_FALSE, flg1, flg2;
34c4762a1bSJed Brown PetscReal ftime;
35c4762a1bSJed Brown PetscInt steps;
36c4762a1bSJed Brown PetscScalar *x_ptr, *y_ptr, *s_ptr;
37c4762a1bSJed Brown Vec lambda[1], q, mu[1];
38c4762a1bSJed Brown PetscInt direction[2];
39c4762a1bSJed Brown PetscBool terminate[2];
40c4762a1bSJed Brown Mat qgrad;
41c4762a1bSJed Brown Mat sp; /* Forward sensitivity matrix */
42c4762a1bSJed Brown SAMethod sa;
43c4762a1bSJed Brown
44c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
45c4762a1bSJed Brown Initialize program
46c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
47327415f7SBarry Smith PetscFunctionBeginUser;
48*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help));
499566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
503c633725SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs");
51c4762a1bSJed Brown
52c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
53c4762a1bSJed Brown Create necessary matrix and vectors
54c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
559566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &ctx.Jac));
569566063dSJacob Faibussowitsch PetscCall(MatSetSizes(ctx.Jac, n, n, PETSC_DETERMINE, PETSC_DETERMINE));
579566063dSJacob Faibussowitsch PetscCall(MatSetType(ctx.Jac, MATDENSE));
589566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(ctx.Jac));
599566063dSJacob Faibussowitsch PetscCall(MatSetUp(ctx.Jac));
609566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(ctx.Jac, &U, NULL));
619566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &ctx.Jacp));
629566063dSJacob Faibussowitsch PetscCall(MatSetSizes(ctx.Jacp, PETSC_DECIDE, PETSC_DECIDE, 2, 1));
639566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(ctx.Jacp));
649566063dSJacob Faibussowitsch PetscCall(MatSetUp(ctx.Jacp));
659566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, 1, 1, NULL, &ctx.DRDP));
669566063dSJacob Faibussowitsch PetscCall(MatSetUp(ctx.DRDP));
679566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, 2, 1, NULL, &ctx.DRDU));
689566063dSJacob Faibussowitsch PetscCall(MatSetUp(ctx.DRDU));
69c4762a1bSJed Brown
70c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
71c4762a1bSJed Brown Set runtime options
72c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
73d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Swing equation options", "");
74c4762a1bSJed Brown {
75c4762a1bSJed Brown ctx.beta = 2;
76c4762a1bSJed Brown ctx.c = 10000.0;
77c4762a1bSJed Brown ctx.u_s = 1.0;
78c4762a1bSJed Brown ctx.omega_s = 1.0;
79c4762a1bSJed Brown ctx.omega_b = 120.0 * PETSC_PI;
80c4762a1bSJed Brown ctx.H = 5.0;
819566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-Inertia", "", "", ctx.H, &ctx.H, NULL));
82c4762a1bSJed Brown ctx.D = 5.0;
839566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-D", "", "", ctx.D, &ctx.D, NULL));
84c4762a1bSJed Brown ctx.E = 1.1378;
85c4762a1bSJed Brown ctx.V = 1.0;
86c4762a1bSJed Brown ctx.X = 0.545;
87c4762a1bSJed Brown ctx.Pmax = ctx.E * ctx.V / ctx.X;
88c4762a1bSJed Brown ctx.Pmax_ini = ctx.Pmax;
899566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-Pmax", "", "", ctx.Pmax, &ctx.Pmax, NULL));
90c4762a1bSJed Brown ctx.Pm = 1.1;
919566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-Pm", "", "", ctx.Pm, &ctx.Pm, NULL));
92c4762a1bSJed Brown ctx.tf = 0.1;
93c4762a1bSJed Brown ctx.tcl = 0.2;
949566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tf", "Time to start fault", "", ctx.tf, &ctx.tf, NULL));
959566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tcl", "Time to end fault", "", ctx.tcl, &ctx.tcl, NULL));
969566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ensemble", "Run ensemble of different initial conditions", "", ensemble, &ensemble, NULL));
97c4762a1bSJed Brown if (ensemble) {
98c4762a1bSJed Brown ctx.tf = -1;
99c4762a1bSJed Brown ctx.tcl = -1;
100c4762a1bSJed Brown }
101c4762a1bSJed Brown
1029566063dSJacob Faibussowitsch PetscCall(VecGetArray(U, &u));
103c4762a1bSJed Brown u[0] = PetscAsinScalar(ctx.Pm / ctx.Pmax);
104c4762a1bSJed Brown u[1] = 1.0;
1059566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-u", "Initial solution", "", u, &n, &flg1));
106c4762a1bSJed Brown n = 2;
1079566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-du", "Perturbation in initial solution", "", du, &n, &flg2));
108c4762a1bSJed Brown u[0] += du[0];
109c4762a1bSJed Brown u[1] += du[1];
1109566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U, &u));
111c4762a1bSJed Brown if (flg1 || flg2) {
112c4762a1bSJed Brown ctx.tf = -1;
113c4762a1bSJed Brown ctx.tcl = -1;
114c4762a1bSJed Brown }
115c4762a1bSJed Brown sa = SA_ADJ;
1169566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-sa_method", "Sensitivity analysis method (adj or tlm)", "", SAMethods, (PetscEnum)sa, (PetscEnum *)&sa, NULL));
117c4762a1bSJed Brown }
118d0609cedSBarry Smith PetscOptionsEnd();
119c4762a1bSJed Brown
120c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
121c4762a1bSJed Brown Create timestepping solver context
122c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1239566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
1249566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
1259566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSBEULER));
1268434afd1SBarry Smith PetscCall(TSSetRHSFunction(ts, NULL, (TSRHSFunctionFn *)RHSFunction, &ctx));
1278434afd1SBarry Smith PetscCall(TSSetRHSJacobian(ts, ctx.Jac, ctx.Jac, (TSRHSJacobianFn *)RHSJacobian, &ctx));
128c4762a1bSJed Brown
129c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
130c4762a1bSJed Brown Set initial conditions
131c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1329566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, U));
133c4762a1bSJed Brown
134c4762a1bSJed Brown /* Set RHS JacobianP */
1359566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobianP(ts, ctx.Jacp, RHSJacobianP, &ctx));
136c4762a1bSJed Brown
1379566063dSJacob Faibussowitsch PetscCall(TSCreateQuadratureTS(ts, PETSC_FALSE, &quadts));
1388434afd1SBarry Smith PetscCall(TSSetRHSFunction(quadts, NULL, (TSRHSFunctionFn *)CostIntegrand, &ctx));
1398434afd1SBarry Smith PetscCall(TSSetRHSJacobian(quadts, ctx.DRDU, ctx.DRDU, (TSRHSJacobianFn *)DRDUJacobianTranspose, &ctx));
1409566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobianP(quadts, ctx.DRDP, DRDPJacobianTranspose, &ctx));
141c4762a1bSJed Brown if (sa == SA_ADJ) {
142c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
143c4762a1bSJed Brown Save trajectory of solution so that TSAdjointSolve() may be used
144c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1459566063dSJacob Faibussowitsch PetscCall(TSSetSaveTrajectory(ts));
1469566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(ctx.Jac, &lambda[0], NULL));
1479566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(ctx.Jacp, &mu[0], NULL));
1489566063dSJacob Faibussowitsch PetscCall(TSSetCostGradients(ts, 1, lambda, mu));
149c4762a1bSJed Brown }
150c4762a1bSJed Brown
151c4762a1bSJed Brown if (sa == SA_TLM) {
152c4762a1bSJed Brown PetscScalar val[2];
153c4762a1bSJed Brown PetscInt row[] = {0, 1}, col[] = {0};
154c4762a1bSJed Brown
1559566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, 1, 1, NULL, &qgrad));
1569566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, 2, 1, NULL, &sp));
1579566063dSJacob Faibussowitsch PetscCall(TSForwardSetSensitivities(ts, 1, sp));
1589566063dSJacob Faibussowitsch PetscCall(TSForwardSetSensitivities(quadts, 1, qgrad));
159c4762a1bSJed Brown val[0] = 1. / PetscSqrtScalar(1. - (ctx.Pm / ctx.Pmax) * (ctx.Pm / ctx.Pmax)) / ctx.Pmax;
160c4762a1bSJed Brown val[1] = 0.0;
1619566063dSJacob Faibussowitsch PetscCall(MatSetValues(sp, 2, row, 1, col, val, INSERT_VALUES));
1629566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(sp, MAT_FINAL_ASSEMBLY));
1639566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(sp, MAT_FINAL_ASSEMBLY));
164c4762a1bSJed Brown }
165c4762a1bSJed Brown
166c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
167c4762a1bSJed Brown Set solver options
168c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1699566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, 1.0));
1709566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
1719566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, 0.03125));
1729566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts));
173c4762a1bSJed Brown
174c4762a1bSJed Brown direction[0] = direction[1] = 1;
175c4762a1bSJed Brown terminate[0] = terminate[1] = PETSC_FALSE;
176c4762a1bSJed Brown
1779566063dSJacob Faibussowitsch PetscCall(TSSetEventHandler(ts, 2, direction, terminate, EventFunction, PostEventFunction, (void *)&ctx));
178c4762a1bSJed Brown
179c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
180c4762a1bSJed Brown Solve nonlinear system
181c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
182c4762a1bSJed Brown if (ensemble) {
183c4762a1bSJed Brown for (du[1] = -2.5; du[1] <= .01; du[1] += .1) {
1849566063dSJacob Faibussowitsch PetscCall(VecGetArray(U, &u));
185c4762a1bSJed Brown u[0] = PetscAsinScalar(ctx.Pm / ctx.Pmax);
186c4762a1bSJed Brown u[1] = ctx.omega_s;
187c4762a1bSJed Brown u[0] += du[0];
188c4762a1bSJed Brown u[1] += du[1];
1899566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U, &u));
1909566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, 0.03125));
1919566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, U));
192c4762a1bSJed Brown }
193c4762a1bSJed Brown } else {
1949566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, U));
195c4762a1bSJed Brown }
1969566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts, &ftime));
1979566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &steps));
198c4762a1bSJed Brown
199c4762a1bSJed Brown if (sa == SA_ADJ) {
200c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
201c4762a1bSJed Brown Adjoint model starts here
202c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
203c4762a1bSJed Brown /* Set initial conditions for the adjoint integration */
2049566063dSJacob Faibussowitsch PetscCall(VecGetArray(lambda[0], &y_ptr));
2059371c9d4SSatish Balay y_ptr[0] = 0.0;
2069371c9d4SSatish Balay y_ptr[1] = 0.0;
2079566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(lambda[0], &y_ptr));
208c4762a1bSJed Brown
2099566063dSJacob Faibussowitsch PetscCall(VecGetArray(mu[0], &x_ptr));
210c4762a1bSJed Brown x_ptr[0] = 0.0;
2119566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(mu[0], &x_ptr));
212c4762a1bSJed Brown
2139566063dSJacob Faibussowitsch PetscCall(TSAdjointSolve(ts));
214c4762a1bSJed Brown
2159566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n lambda: d[Psi(tf)]/d[phi0] d[Psi(tf)]/d[omega0]\n"));
2169566063dSJacob Faibussowitsch PetscCall(VecView(lambda[0], PETSC_VIEWER_STDOUT_WORLD));
2179566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n mu: d[Psi(tf)]/d[pm]\n"));
2189566063dSJacob Faibussowitsch PetscCall(VecView(mu[0], PETSC_VIEWER_STDOUT_WORLD));
2199566063dSJacob Faibussowitsch PetscCall(TSGetCostIntegral(ts, &q));
2209566063dSJacob Faibussowitsch PetscCall(VecGetArray(q, &x_ptr));
2219566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n cost function=%g\n", (double)(x_ptr[0] - ctx.Pm)));
2229566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(q, &x_ptr));
2239566063dSJacob Faibussowitsch PetscCall(ComputeSensiP(lambda[0], mu[0], &ctx));
2249566063dSJacob Faibussowitsch PetscCall(VecGetArray(mu[0], &x_ptr));
2259566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n gradient=%g\n", (double)x_ptr[0]));
2269566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(mu[0], &x_ptr));
2279566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lambda[0]));
2289566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mu[0]));
229c4762a1bSJed Brown }
230c4762a1bSJed Brown if (sa == SA_TLM) {
2319566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n trajectory sensitivity: d[phi(tf)]/d[pm] d[omega(tf)]/d[pm]\n"));
2329566063dSJacob Faibussowitsch PetscCall(MatView(sp, PETSC_VIEWER_STDOUT_WORLD));
2339566063dSJacob Faibussowitsch PetscCall(TSGetCostIntegral(ts, &q));
2349566063dSJacob Faibussowitsch PetscCall(VecGetArray(q, &s_ptr));
2359566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n cost function=%g\n", (double)(s_ptr[0] - ctx.Pm)));
2369566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(q, &s_ptr));
2379566063dSJacob Faibussowitsch PetscCall(MatDenseGetArray(qgrad, &s_ptr));
2389566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n gradient=%g\n", (double)s_ptr[0]));
2399566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArray(qgrad, &s_ptr));
2409566063dSJacob Faibussowitsch PetscCall(MatDestroy(&qgrad));
2419566063dSJacob Faibussowitsch PetscCall(MatDestroy(&sp));
242c4762a1bSJed Brown }
243c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
244c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they are no longer needed.
245c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2469566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ctx.Jac));
2479566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ctx.Jacp));
2489566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ctx.DRDU));
2499566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ctx.DRDP));
2509566063dSJacob Faibussowitsch PetscCall(VecDestroy(&U));
2519566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts));
2529566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
253b122ec5aSJacob Faibussowitsch return 0;
254c4762a1bSJed Brown }
255c4762a1bSJed Brown
256c4762a1bSJed Brown /*TEST
257c4762a1bSJed Brown
258c4762a1bSJed Brown build:
259c4762a1bSJed Brown requires: !complex !single
260c4762a1bSJed Brown
261c4762a1bSJed Brown test:
262c4762a1bSJed Brown args: -sa_method adj -viewer_binary_skip_info -ts_type cn -pc_type lu
263c4762a1bSJed Brown
264c4762a1bSJed Brown test:
265c4762a1bSJed Brown suffix: 2
266c4762a1bSJed Brown args: -sa_method tlm -ts_type cn -pc_type lu
267c4762a1bSJed Brown
268c4762a1bSJed Brown test:
269c4762a1bSJed Brown suffix: 3
270c4762a1bSJed Brown args: -sa_method adj -ts_type rk -ts_rk_type 2a -ts_adapt_type dsp
271c4762a1bSJed Brown
272c4762a1bSJed Brown TEST*/
273