xref: /petsc/src/ts/tutorials/power_grid/ex3.h (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1*9371c9d4SSatish Balay typedef enum {
2*9371c9d4SSatish Balay   SA_ADJ,
3*9371c9d4SSatish Balay   SA_TLM
4*9371c9d4SSatish Balay } SAMethod;
5c4762a1bSJed Brown static const char *const SAMethods[] = {"ADJ", "TLM", "SAMethod", "SA_", 0};
6c4762a1bSJed Brown 
7c4762a1bSJed Brown typedef struct {
8c4762a1bSJed Brown   PetscScalar H, D, omega_b, omega_s, Pmax, Pmax_ini, Pm, E, V, X, u_s, c;
9c4762a1bSJed Brown   PetscInt    beta;
10c4762a1bSJed Brown   PetscReal   tf, tcl;
11c4762a1bSJed Brown   /* Solver context */
12c4762a1bSJed Brown   TS          ts, quadts;
13c4762a1bSJed Brown   Vec         U;    /* solution will be stored here */
14c4762a1bSJed Brown   Mat         Jac;  /* Jacobian matrix */
15c4762a1bSJed Brown   Mat         Jacp; /* Jacobianp matrix */
16c4762a1bSJed Brown   Mat         DRDU, DRDP;
17c4762a1bSJed Brown   SAMethod    sa;
18c4762a1bSJed Brown } AppCtx;
19c4762a1bSJed Brown 
20c4762a1bSJed Brown /* Event check */
21*9371c9d4SSatish Balay PetscErrorCode EventFunction(TS ts, PetscReal t, Vec X, PetscScalar *fvalue, void *ctx) {
22c4762a1bSJed Brown   AppCtx *user = (AppCtx *)ctx;
23c4762a1bSJed Brown 
24c4762a1bSJed Brown   PetscFunctionBegin;
25c4762a1bSJed Brown   /* Event for fault-on time */
26c4762a1bSJed Brown   fvalue[0] = t - user->tf;
27c4762a1bSJed Brown   /* Event for fault-off time */
28c4762a1bSJed Brown   fvalue[1] = t - user->tcl;
29c4762a1bSJed Brown 
30c4762a1bSJed Brown   PetscFunctionReturn(0);
31c4762a1bSJed Brown }
32c4762a1bSJed Brown 
33*9371c9d4SSatish Balay PetscErrorCode PostEventFunction(TS ts, PetscInt nevents, PetscInt event_list[], PetscReal t, Vec X, PetscBool forwardsolve, void *ctx) {
34c4762a1bSJed Brown   AppCtx *user = (AppCtx *)ctx;
35c4762a1bSJed Brown 
36c4762a1bSJed Brown   PetscFunctionBegin;
37c4762a1bSJed Brown   if (event_list[0] == 0) {
38c4762a1bSJed Brown     if (forwardsolve) user->Pmax = 0.0; /* Apply disturbance - this is done by setting Pmax = 0 */
39c4762a1bSJed Brown     else user->Pmax = user->Pmax_ini;   /* Going backward, reversal of event */
40c4762a1bSJed Brown   } else if (event_list[0] == 1) {
41c4762a1bSJed Brown     if (forwardsolve) user->Pmax = user->Pmax_ini; /* Remove the fault  - this is done by setting Pmax = Pmax_ini */
42c4762a1bSJed Brown     else user->Pmax = 0.0;                         /* Going backward, reversal of event */
43c4762a1bSJed Brown   }
449566063dSJacob Faibussowitsch   PetscCall(TSRestartStep(ts)); /* Must set restart flag to ture, otherwise methods with FSAL will fail */
45c4762a1bSJed Brown   PetscFunctionReturn(0);
46c4762a1bSJed Brown }
47c4762a1bSJed Brown 
48c4762a1bSJed Brown /*
49c4762a1bSJed Brown      Defines the ODE passed to the ODE solver
50c4762a1bSJed Brown */
51*9371c9d4SSatish Balay PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec U, Vec F, AppCtx *ctx) {
52c4762a1bSJed Brown   PetscScalar       *f, Pmax;
53c4762a1bSJed Brown   const PetscScalar *u;
54c4762a1bSJed Brown 
55c4762a1bSJed Brown   PetscFunctionBegin;
56c4762a1bSJed Brown   /*  The next three lines allow us to access the entries of the vectors directly */
579566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
589566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
59c4762a1bSJed Brown   Pmax = ctx->Pmax;
60c4762a1bSJed Brown   f[0] = ctx->omega_b * (u[1] - ctx->omega_s);
61c4762a1bSJed Brown   f[1] = ctx->omega_s / (2.0 * ctx->H) * (ctx->Pm - Pmax * PetscSinScalar(u[0]) - ctx->D * (u[1] - ctx->omega_s));
62c4762a1bSJed Brown 
639566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
649566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
65c4762a1bSJed Brown   PetscFunctionReturn(0);
66c4762a1bSJed Brown }
67c4762a1bSJed Brown 
68c4762a1bSJed Brown /*
69c4762a1bSJed Brown      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetRHSJacobian() for the meaning of a and the Jacobian.
70c4762a1bSJed Brown */
71*9371c9d4SSatish Balay PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec U, Mat A, Mat B, AppCtx *ctx) {
72c4762a1bSJed Brown   PetscInt           rowcol[] = {0, 1};
73c4762a1bSJed Brown   PetscScalar        J[2][2], Pmax;
74c4762a1bSJed Brown   const PetscScalar *u;
75c4762a1bSJed Brown 
76c4762a1bSJed Brown   PetscFunctionBegin;
779566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
78c4762a1bSJed Brown   Pmax    = ctx->Pmax;
79c4762a1bSJed Brown   J[0][0] = 0;
80c4762a1bSJed Brown   J[0][1] = ctx->omega_b;
81c4762a1bSJed Brown   J[1][0] = -ctx->omega_s / (2.0 * ctx->H) * Pmax * PetscCosScalar(u[0]);
82c4762a1bSJed Brown   J[1][1] = -ctx->omega_s / (2.0 * ctx->H) * ctx->D;
839566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES));
849566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
85c4762a1bSJed Brown 
869566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
879566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
88c4762a1bSJed Brown   if (A != B) {
899566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
909566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
91c4762a1bSJed Brown   }
92c4762a1bSJed Brown   PetscFunctionReturn(0);
93c4762a1bSJed Brown }
94c4762a1bSJed Brown 
95c4762a1bSJed Brown /*
96c4762a1bSJed Brown      Defines the ODE passed to the ODE solver
97c4762a1bSJed Brown */
98*9371c9d4SSatish Balay PetscErrorCode IFunction(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, AppCtx *ctx) {
99c4762a1bSJed Brown   PetscScalar       *f, Pmax;
100c4762a1bSJed Brown   const PetscScalar *u, *udot;
101c4762a1bSJed Brown 
102c4762a1bSJed Brown   PetscFunctionBegin;
103c4762a1bSJed Brown   /*  The next three lines allow us to access the entries of the vectors directly */
1049566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1059566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Udot, &udot));
1069566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
107c4762a1bSJed Brown   Pmax = ctx->Pmax;
108c4762a1bSJed Brown   f[0] = udot[0] - ctx->omega_b * (u[1] - ctx->omega_s);
109c4762a1bSJed Brown   f[1] = 2.0 * ctx->H / ctx->omega_s * udot[1] + Pmax * PetscSinScalar(u[0]) + ctx->D * (u[1] - ctx->omega_s) - ctx->Pm;
110c4762a1bSJed Brown 
1119566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1129566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Udot, &udot));
1139566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
114c4762a1bSJed Brown   PetscFunctionReturn(0);
115c4762a1bSJed Brown }
116c4762a1bSJed Brown 
117c4762a1bSJed Brown /*
118c4762a1bSJed Brown      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
119c4762a1bSJed Brown */
120*9371c9d4SSatish Balay PetscErrorCode IJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat B, AppCtx *ctx) {
121c4762a1bSJed Brown   PetscInt           rowcol[] = {0, 1};
122c4762a1bSJed Brown   PetscScalar        J[2][2], Pmax;
123c4762a1bSJed Brown   const PetscScalar *u, *udot;
124c4762a1bSJed Brown 
125c4762a1bSJed Brown   PetscFunctionBegin;
1269566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1279566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Udot, &udot));
128c4762a1bSJed Brown   Pmax    = ctx->Pmax;
129*9371c9d4SSatish Balay   J[0][0] = a;
130*9371c9d4SSatish Balay   J[0][1] = -ctx->omega_b;
131*9371c9d4SSatish Balay   J[1][1] = 2.0 * ctx->H / ctx->omega_s * a + ctx->D;
132*9371c9d4SSatish Balay   J[1][0] = Pmax * PetscCosScalar(u[0]);
133c4762a1bSJed Brown 
1349566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES));
1359566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1369566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Udot, &udot));
137c4762a1bSJed Brown 
1389566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1399566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
140c4762a1bSJed Brown   if (A != B) {
1419566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
1429566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
143c4762a1bSJed Brown   }
144c4762a1bSJed Brown   PetscFunctionReturn(0);
145c4762a1bSJed Brown }
146c4762a1bSJed Brown 
147*9371c9d4SSatish Balay PetscErrorCode RHSJacobianP(TS ts, PetscReal t, Vec X, Mat A, void *ctx0) {
148c4762a1bSJed Brown   PetscInt     row[] = {0, 1}, col[] = {0};
149c4762a1bSJed Brown   PetscScalar *x, J[2][1];
150c4762a1bSJed Brown   AppCtx      *ctx = (AppCtx *)ctx0;
151c4762a1bSJed Brown 
152c4762a1bSJed Brown   PetscFunctionBeginUser;
1539566063dSJacob Faibussowitsch   PetscCall(VecGetArray(X, &x));
154c4762a1bSJed Brown   J[0][0] = 0;
155c4762a1bSJed Brown   J[1][0] = ctx->omega_s / (2.0 * ctx->H);
1569566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 2, row, 1, col, &J[0][0], INSERT_VALUES));
157c4762a1bSJed Brown 
1589566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1599566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
160c4762a1bSJed Brown   PetscFunctionReturn(0);
161c4762a1bSJed Brown }
162c4762a1bSJed Brown 
163*9371c9d4SSatish Balay PetscErrorCode CostIntegrand(TS ts, PetscReal t, Vec U, Vec R, AppCtx *ctx) {
164c4762a1bSJed Brown   PetscScalar       *r;
165c4762a1bSJed Brown   const PetscScalar *u;
166c4762a1bSJed Brown 
167c4762a1bSJed Brown   PetscFunctionBegin;
1689566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1699566063dSJacob Faibussowitsch   PetscCall(VecGetArray(R, &r));
1702f613bf5SBarry Smith   r[0] = ctx->c * PetscPowScalarInt(PetscMax(0., u[0] - ctx->u_s), ctx->beta);
1719566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(R, &r));
1729566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
173c4762a1bSJed Brown   PetscFunctionReturn(0);
174c4762a1bSJed Brown }
175c4762a1bSJed Brown 
176c4762a1bSJed Brown /* Transpose of DRDU */
177*9371c9d4SSatish Balay PetscErrorCode DRDUJacobianTranspose(TS ts, PetscReal t, Vec U, Mat DRDU, Mat B, AppCtx *ctx) {
178c4762a1bSJed Brown   PetscScalar        ru[2];
179c4762a1bSJed Brown   PetscInt           row[] = {0, 1}, col[] = {0};
180c4762a1bSJed Brown   const PetscScalar *u;
181c4762a1bSJed Brown 
182c4762a1bSJed Brown   PetscFunctionBegin;
1839566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1842f613bf5SBarry Smith   ru[0] = ctx->c * ctx->beta * PetscPowScalarInt(PetscMax(0., u[0] - ctx->u_s), ctx->beta - 1);
185c4762a1bSJed Brown   ru[1] = 0;
1869566063dSJacob Faibussowitsch   PetscCall(MatSetValues(DRDU, 2, row, 1, col, ru, INSERT_VALUES));
1879566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1889566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(DRDU, MAT_FINAL_ASSEMBLY));
1899566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(DRDU, MAT_FINAL_ASSEMBLY));
190c4762a1bSJed Brown   PetscFunctionReturn(0);
191c4762a1bSJed Brown }
192c4762a1bSJed Brown 
193*9371c9d4SSatish Balay PetscErrorCode DRDPJacobianTranspose(TS ts, PetscReal t, Vec U, Mat DRDP, void *ctx) {
194c4762a1bSJed Brown   PetscFunctionBegin;
1959566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(DRDP));
1969566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(DRDP, MAT_FINAL_ASSEMBLY));
1979566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(DRDP, MAT_FINAL_ASSEMBLY));
198c4762a1bSJed Brown   PetscFunctionReturn(0);
199c4762a1bSJed Brown }
200c4762a1bSJed Brown 
201*9371c9d4SSatish Balay PetscErrorCode ComputeSensiP(Vec lambda, Vec mu, AppCtx *ctx) {
202c4762a1bSJed Brown   PetscScalar       *y, sensip;
203c4762a1bSJed Brown   const PetscScalar *x;
204c4762a1bSJed Brown 
205c4762a1bSJed Brown   PetscFunctionBegin;
2069566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(lambda, &x));
2079566063dSJacob Faibussowitsch   PetscCall(VecGetArray(mu, &y));
208c4762a1bSJed Brown   sensip = 1. / PetscSqrtScalar(1. - (ctx->Pm / ctx->Pmax) * (ctx->Pm / ctx->Pmax)) / ctx->Pmax * x[0] + y[0];
209c4762a1bSJed Brown   y[0]   = sensip;
2109566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(mu, &y));
2119566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(lambda, &x));
212c4762a1bSJed Brown   PetscFunctionReturn(0);
213c4762a1bSJed Brown }
214