xref: /petsc/src/ts/tutorials/power_grid/ex3.h (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
19371c9d4SSatish Balay typedef enum {
29371c9d4SSatish Balay   SA_ADJ,
39371c9d4SSatish Balay   SA_TLM
49371c9d4SSatish 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*d71ae5a4SJacob Faibussowitsch PetscErrorCode EventFunction(TS ts, PetscReal t, Vec X, PetscScalar *fvalue, void *ctx)
22*d71ae5a4SJacob Faibussowitsch {
23c4762a1bSJed Brown   AppCtx *user = (AppCtx *)ctx;
24c4762a1bSJed Brown 
25c4762a1bSJed Brown   PetscFunctionBegin;
26c4762a1bSJed Brown   /* Event for fault-on time */
27c4762a1bSJed Brown   fvalue[0] = t - user->tf;
28c4762a1bSJed Brown   /* Event for fault-off time */
29c4762a1bSJed Brown   fvalue[1] = t - user->tcl;
30c4762a1bSJed Brown 
31c4762a1bSJed Brown   PetscFunctionReturn(0);
32c4762a1bSJed Brown }
33c4762a1bSJed Brown 
34*d71ae5a4SJacob Faibussowitsch PetscErrorCode PostEventFunction(TS ts, PetscInt nevents, PetscInt event_list[], PetscReal t, Vec X, PetscBool forwardsolve, void *ctx)
35*d71ae5a4SJacob Faibussowitsch {
36c4762a1bSJed Brown   AppCtx *user = (AppCtx *)ctx;
37c4762a1bSJed Brown 
38c4762a1bSJed Brown   PetscFunctionBegin;
39c4762a1bSJed Brown   if (event_list[0] == 0) {
40c4762a1bSJed Brown     if (forwardsolve) user->Pmax = 0.0; /* Apply disturbance - this is done by setting Pmax = 0 */
41c4762a1bSJed Brown     else user->Pmax = user->Pmax_ini;   /* Going backward, reversal of event */
42c4762a1bSJed Brown   } else if (event_list[0] == 1) {
43c4762a1bSJed Brown     if (forwardsolve) user->Pmax = user->Pmax_ini; /* Remove the fault  - this is done by setting Pmax = Pmax_ini */
44c4762a1bSJed Brown     else user->Pmax = 0.0;                         /* Going backward, reversal of event */
45c4762a1bSJed Brown   }
469566063dSJacob Faibussowitsch   PetscCall(TSRestartStep(ts)); /* Must set restart flag to ture, otherwise methods with FSAL will fail */
47c4762a1bSJed Brown   PetscFunctionReturn(0);
48c4762a1bSJed Brown }
49c4762a1bSJed Brown 
50c4762a1bSJed Brown /*
51c4762a1bSJed Brown      Defines the ODE passed to the ODE solver
52c4762a1bSJed Brown */
53*d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec U, Vec F, AppCtx *ctx)
54*d71ae5a4SJacob Faibussowitsch {
55c4762a1bSJed Brown   PetscScalar       *f, Pmax;
56c4762a1bSJed Brown   const PetscScalar *u;
57c4762a1bSJed Brown 
58c4762a1bSJed Brown   PetscFunctionBegin;
59c4762a1bSJed Brown   /*  The next three lines allow us to access the entries of the vectors directly */
609566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
619566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
62c4762a1bSJed Brown   Pmax = ctx->Pmax;
63c4762a1bSJed Brown   f[0] = ctx->omega_b * (u[1] - ctx->omega_s);
64c4762a1bSJed Brown   f[1] = ctx->omega_s / (2.0 * ctx->H) * (ctx->Pm - Pmax * PetscSinScalar(u[0]) - ctx->D * (u[1] - ctx->omega_s));
65c4762a1bSJed Brown 
669566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
679566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
68c4762a1bSJed Brown   PetscFunctionReturn(0);
69c4762a1bSJed Brown }
70c4762a1bSJed Brown 
71c4762a1bSJed Brown /*
72c4762a1bSJed Brown      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetRHSJacobian() for the meaning of a and the Jacobian.
73c4762a1bSJed Brown */
74*d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec U, Mat A, Mat B, AppCtx *ctx)
75*d71ae5a4SJacob Faibussowitsch {
76c4762a1bSJed Brown   PetscInt           rowcol[] = {0, 1};
77c4762a1bSJed Brown   PetscScalar        J[2][2], Pmax;
78c4762a1bSJed Brown   const PetscScalar *u;
79c4762a1bSJed Brown 
80c4762a1bSJed Brown   PetscFunctionBegin;
819566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
82c4762a1bSJed Brown   Pmax    = ctx->Pmax;
83c4762a1bSJed Brown   J[0][0] = 0;
84c4762a1bSJed Brown   J[0][1] = ctx->omega_b;
85c4762a1bSJed Brown   J[1][0] = -ctx->omega_s / (2.0 * ctx->H) * Pmax * PetscCosScalar(u[0]);
86c4762a1bSJed Brown   J[1][1] = -ctx->omega_s / (2.0 * ctx->H) * ctx->D;
879566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES));
889566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
89c4762a1bSJed Brown 
909566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
919566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
92c4762a1bSJed Brown   if (A != B) {
939566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
949566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
95c4762a1bSJed Brown   }
96c4762a1bSJed Brown   PetscFunctionReturn(0);
97c4762a1bSJed Brown }
98c4762a1bSJed Brown 
99c4762a1bSJed Brown /*
100c4762a1bSJed Brown      Defines the ODE passed to the ODE solver
101c4762a1bSJed Brown */
102*d71ae5a4SJacob Faibussowitsch PetscErrorCode IFunction(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, AppCtx *ctx)
103*d71ae5a4SJacob Faibussowitsch {
104c4762a1bSJed Brown   PetscScalar       *f, Pmax;
105c4762a1bSJed Brown   const PetscScalar *u, *udot;
106c4762a1bSJed Brown 
107c4762a1bSJed Brown   PetscFunctionBegin;
108c4762a1bSJed Brown   /*  The next three lines allow us to access the entries of the vectors directly */
1099566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1109566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Udot, &udot));
1119566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
112c4762a1bSJed Brown   Pmax = ctx->Pmax;
113c4762a1bSJed Brown   f[0] = udot[0] - ctx->omega_b * (u[1] - ctx->omega_s);
114c4762a1bSJed 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;
115c4762a1bSJed Brown 
1169566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1179566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Udot, &udot));
1189566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
119c4762a1bSJed Brown   PetscFunctionReturn(0);
120c4762a1bSJed Brown }
121c4762a1bSJed Brown 
122c4762a1bSJed Brown /*
123c4762a1bSJed Brown      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
124c4762a1bSJed Brown */
125*d71ae5a4SJacob Faibussowitsch PetscErrorCode IJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat B, AppCtx *ctx)
126*d71ae5a4SJacob Faibussowitsch {
127c4762a1bSJed Brown   PetscInt           rowcol[] = {0, 1};
128c4762a1bSJed Brown   PetscScalar        J[2][2], Pmax;
129c4762a1bSJed Brown   const PetscScalar *u, *udot;
130c4762a1bSJed Brown 
131c4762a1bSJed Brown   PetscFunctionBegin;
1329566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1339566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Udot, &udot));
134c4762a1bSJed Brown   Pmax    = ctx->Pmax;
1359371c9d4SSatish Balay   J[0][0] = a;
1369371c9d4SSatish Balay   J[0][1] = -ctx->omega_b;
1379371c9d4SSatish Balay   J[1][1] = 2.0 * ctx->H / ctx->omega_s * a + ctx->D;
1389371c9d4SSatish Balay   J[1][0] = Pmax * PetscCosScalar(u[0]);
139c4762a1bSJed Brown 
1409566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES));
1419566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1429566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Udot, &udot));
143c4762a1bSJed Brown 
1449566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1459566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
146c4762a1bSJed Brown   if (A != B) {
1479566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
1489566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
149c4762a1bSJed Brown   }
150c4762a1bSJed Brown   PetscFunctionReturn(0);
151c4762a1bSJed Brown }
152c4762a1bSJed Brown 
153*d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSJacobianP(TS ts, PetscReal t, Vec X, Mat A, void *ctx0)
154*d71ae5a4SJacob Faibussowitsch {
155c4762a1bSJed Brown   PetscInt     row[] = {0, 1}, col[] = {0};
156c4762a1bSJed Brown   PetscScalar *x, J[2][1];
157c4762a1bSJed Brown   AppCtx      *ctx = (AppCtx *)ctx0;
158c4762a1bSJed Brown 
159c4762a1bSJed Brown   PetscFunctionBeginUser;
1609566063dSJacob Faibussowitsch   PetscCall(VecGetArray(X, &x));
161c4762a1bSJed Brown   J[0][0] = 0;
162c4762a1bSJed Brown   J[1][0] = ctx->omega_s / (2.0 * ctx->H);
1639566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 2, row, 1, col, &J[0][0], INSERT_VALUES));
164c4762a1bSJed Brown 
1659566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1669566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
167c4762a1bSJed Brown   PetscFunctionReturn(0);
168c4762a1bSJed Brown }
169c4762a1bSJed Brown 
170*d71ae5a4SJacob Faibussowitsch PetscErrorCode CostIntegrand(TS ts, PetscReal t, Vec U, Vec R, AppCtx *ctx)
171*d71ae5a4SJacob Faibussowitsch {
172c4762a1bSJed Brown   PetscScalar       *r;
173c4762a1bSJed Brown   const PetscScalar *u;
174c4762a1bSJed Brown 
175c4762a1bSJed Brown   PetscFunctionBegin;
1769566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1779566063dSJacob Faibussowitsch   PetscCall(VecGetArray(R, &r));
1782f613bf5SBarry Smith   r[0] = ctx->c * PetscPowScalarInt(PetscMax(0., u[0] - ctx->u_s), ctx->beta);
1799566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(R, &r));
1809566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
181c4762a1bSJed Brown   PetscFunctionReturn(0);
182c4762a1bSJed Brown }
183c4762a1bSJed Brown 
184c4762a1bSJed Brown /* Transpose of DRDU */
185*d71ae5a4SJacob Faibussowitsch PetscErrorCode DRDUJacobianTranspose(TS ts, PetscReal t, Vec U, Mat DRDU, Mat B, AppCtx *ctx)
186*d71ae5a4SJacob Faibussowitsch {
187c4762a1bSJed Brown   PetscScalar        ru[2];
188c4762a1bSJed Brown   PetscInt           row[] = {0, 1}, col[] = {0};
189c4762a1bSJed Brown   const PetscScalar *u;
190c4762a1bSJed Brown 
191c4762a1bSJed Brown   PetscFunctionBegin;
1929566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1932f613bf5SBarry Smith   ru[0] = ctx->c * ctx->beta * PetscPowScalarInt(PetscMax(0., u[0] - ctx->u_s), ctx->beta - 1);
194c4762a1bSJed Brown   ru[1] = 0;
1959566063dSJacob Faibussowitsch   PetscCall(MatSetValues(DRDU, 2, row, 1, col, ru, INSERT_VALUES));
1969566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1979566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(DRDU, MAT_FINAL_ASSEMBLY));
1989566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(DRDU, MAT_FINAL_ASSEMBLY));
199c4762a1bSJed Brown   PetscFunctionReturn(0);
200c4762a1bSJed Brown }
201c4762a1bSJed Brown 
202*d71ae5a4SJacob Faibussowitsch PetscErrorCode DRDPJacobianTranspose(TS ts, PetscReal t, Vec U, Mat DRDP, void *ctx)
203*d71ae5a4SJacob Faibussowitsch {
204c4762a1bSJed Brown   PetscFunctionBegin;
2059566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(DRDP));
2069566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(DRDP, MAT_FINAL_ASSEMBLY));
2079566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(DRDP, MAT_FINAL_ASSEMBLY));
208c4762a1bSJed Brown   PetscFunctionReturn(0);
209c4762a1bSJed Brown }
210c4762a1bSJed Brown 
211*d71ae5a4SJacob Faibussowitsch PetscErrorCode ComputeSensiP(Vec lambda, Vec mu, AppCtx *ctx)
212*d71ae5a4SJacob Faibussowitsch {
213c4762a1bSJed Brown   PetscScalar       *y, sensip;
214c4762a1bSJed Brown   const PetscScalar *x;
215c4762a1bSJed Brown 
216c4762a1bSJed Brown   PetscFunctionBegin;
2179566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(lambda, &x));
2189566063dSJacob Faibussowitsch   PetscCall(VecGetArray(mu, &y));
219c4762a1bSJed Brown   sensip = 1. / PetscSqrtScalar(1. - (ctx->Pm / ctx->Pmax) * (ctx->Pm / ctx->Pmax)) / ctx->Pmax * x[0] + y[0];
220c4762a1bSJed Brown   y[0]   = sensip;
2219566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(mu, &y));
2229566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(lambda, &x));
223c4762a1bSJed Brown   PetscFunctionReturn(0);
224c4762a1bSJed Brown }
225