xref: /petsc/src/ts/tutorials/power_grid/ex3.h (revision a496304597bacff3545e802853d69e8765312868)
1*a4963045SJacob Faibussowitsch #pragma once
23ea99036SJacob Faibussowitsch 
39371c9d4SSatish Balay typedef enum {
49371c9d4SSatish Balay   SA_ADJ,
59371c9d4SSatish Balay   SA_TLM
69371c9d4SSatish Balay } SAMethod;
7c4762a1bSJed Brown static const char *const SAMethods[] = {"ADJ", "TLM", "SAMethod", "SA_", 0};
8c4762a1bSJed Brown 
9c4762a1bSJed Brown typedef struct {
10c4762a1bSJed Brown   PetscScalar H, D, omega_b, omega_s, Pmax, Pmax_ini, Pm, E, V, X, u_s, c;
11c4762a1bSJed Brown   PetscInt    beta;
12c4762a1bSJed Brown   PetscReal   tf, tcl;
13c4762a1bSJed Brown   /* Solver context */
14c4762a1bSJed Brown   TS       ts, quadts;
15c4762a1bSJed Brown   Vec      U;    /* solution will be stored here */
16c4762a1bSJed Brown   Mat      Jac;  /* Jacobian matrix */
17c4762a1bSJed Brown   Mat      Jacp; /* Jacobianp matrix */
18c4762a1bSJed Brown   Mat      DRDU, DRDP;
19c4762a1bSJed Brown   SAMethod sa;
20c4762a1bSJed Brown } AppCtx;
21c4762a1bSJed Brown 
22c4762a1bSJed Brown /* Event check */
23d71ae5a4SJacob Faibussowitsch PetscErrorCode EventFunction(TS ts, PetscReal t, Vec X, PetscScalar *fvalue, void *ctx)
24d71ae5a4SJacob Faibussowitsch {
25c4762a1bSJed Brown   AppCtx *user = (AppCtx *)ctx;
26c4762a1bSJed Brown 
27c4762a1bSJed Brown   PetscFunctionBegin;
28c4762a1bSJed Brown   /* Event for fault-on time */
29c4762a1bSJed Brown   fvalue[0] = t - user->tf;
30c4762a1bSJed Brown   /* Event for fault-off time */
31c4762a1bSJed Brown   fvalue[1] = t - user->tcl;
32c4762a1bSJed Brown 
333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
34c4762a1bSJed Brown }
35c4762a1bSJed Brown 
36d71ae5a4SJacob Faibussowitsch PetscErrorCode PostEventFunction(TS ts, PetscInt nevents, PetscInt event_list[], PetscReal t, Vec X, PetscBool forwardsolve, void *ctx)
37d71ae5a4SJacob Faibussowitsch {
38c4762a1bSJed Brown   AppCtx *user = (AppCtx *)ctx;
39c4762a1bSJed Brown 
40c4762a1bSJed Brown   PetscFunctionBegin;
41c4762a1bSJed Brown   if (event_list[0] == 0) {
42c4762a1bSJed Brown     if (forwardsolve) user->Pmax = 0.0; /* Apply disturbance - this is done by setting Pmax = 0 */
43c4762a1bSJed Brown     else user->Pmax = user->Pmax_ini;   /* Going backward, reversal of event */
44c4762a1bSJed Brown   } else if (event_list[0] == 1) {
45c4762a1bSJed Brown     if (forwardsolve) user->Pmax = user->Pmax_ini; /* Remove the fault  - this is done by setting Pmax = Pmax_ini */
46c4762a1bSJed Brown     else user->Pmax = 0.0;                         /* Going backward, reversal of event */
47c4762a1bSJed Brown   }
4869d47153SPierre Jolivet   PetscCall(TSRestartStep(ts)); /* Must set restart flag to true, otherwise methods with FSAL will fail */
493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
50c4762a1bSJed Brown }
51c4762a1bSJed Brown 
52c4762a1bSJed Brown /*
53c4762a1bSJed Brown      Defines the ODE passed to the ODE solver
54c4762a1bSJed Brown */
55d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec U, Vec F, AppCtx *ctx)
56d71ae5a4SJacob Faibussowitsch {
57c4762a1bSJed Brown   PetscScalar       *f, Pmax;
58c4762a1bSJed Brown   const PetscScalar *u;
59c4762a1bSJed Brown 
60c4762a1bSJed Brown   PetscFunctionBegin;
61c4762a1bSJed Brown   /*  The next three lines allow us to access the entries of the vectors directly */
629566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
639566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
64c4762a1bSJed Brown   Pmax = ctx->Pmax;
65c4762a1bSJed Brown   f[0] = ctx->omega_b * (u[1] - ctx->omega_s);
66c4762a1bSJed Brown   f[1] = ctx->omega_s / (2.0 * ctx->H) * (ctx->Pm - Pmax * PetscSinScalar(u[0]) - ctx->D * (u[1] - ctx->omega_s));
67c4762a1bSJed Brown 
689566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
699566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
71c4762a1bSJed Brown }
72c4762a1bSJed Brown 
73c4762a1bSJed Brown /*
74c4762a1bSJed Brown      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetRHSJacobian() for the meaning of a and the Jacobian.
75c4762a1bSJed Brown */
76d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec U, Mat A, Mat B, AppCtx *ctx)
77d71ae5a4SJacob Faibussowitsch {
78c4762a1bSJed Brown   PetscInt           rowcol[] = {0, 1};
79c4762a1bSJed Brown   PetscScalar        J[2][2], Pmax;
80c4762a1bSJed Brown   const PetscScalar *u;
81c4762a1bSJed Brown 
82c4762a1bSJed Brown   PetscFunctionBegin;
839566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
84c4762a1bSJed Brown   Pmax    = ctx->Pmax;
85c4762a1bSJed Brown   J[0][0] = 0;
86c4762a1bSJed Brown   J[0][1] = ctx->omega_b;
87c4762a1bSJed Brown   J[1][0] = -ctx->omega_s / (2.0 * ctx->H) * Pmax * PetscCosScalar(u[0]);
88c4762a1bSJed Brown   J[1][1] = -ctx->omega_s / (2.0 * ctx->H) * ctx->D;
899566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES));
909566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
91c4762a1bSJed Brown 
929566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
939566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
94c4762a1bSJed Brown   if (A != B) {
959566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
969566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
97c4762a1bSJed Brown   }
983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
99c4762a1bSJed Brown }
100c4762a1bSJed Brown 
101c4762a1bSJed Brown /*
102c4762a1bSJed Brown      Defines the ODE passed to the ODE solver
103c4762a1bSJed Brown */
104d71ae5a4SJacob Faibussowitsch PetscErrorCode IFunction(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, AppCtx *ctx)
105d71ae5a4SJacob Faibussowitsch {
106c4762a1bSJed Brown   PetscScalar       *f, Pmax;
107c4762a1bSJed Brown   const PetscScalar *u, *udot;
108c4762a1bSJed Brown 
109c4762a1bSJed Brown   PetscFunctionBegin;
110c4762a1bSJed Brown   /*  The next three lines allow us to access the entries of the vectors directly */
1119566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1129566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Udot, &udot));
1139566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
114c4762a1bSJed Brown   Pmax = ctx->Pmax;
115c4762a1bSJed Brown   f[0] = udot[0] - ctx->omega_b * (u[1] - ctx->omega_s);
116c4762a1bSJed 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;
117c4762a1bSJed Brown 
1189566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1199566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Udot, &udot));
1209566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
1213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
122c4762a1bSJed Brown }
123c4762a1bSJed Brown 
124c4762a1bSJed Brown /*
125c4762a1bSJed Brown      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
126c4762a1bSJed Brown */
127d71ae5a4SJacob Faibussowitsch PetscErrorCode IJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat B, AppCtx *ctx)
128d71ae5a4SJacob Faibussowitsch {
129c4762a1bSJed Brown   PetscInt           rowcol[] = {0, 1};
130c4762a1bSJed Brown   PetscScalar        J[2][2], Pmax;
131c4762a1bSJed Brown   const PetscScalar *u, *udot;
132c4762a1bSJed Brown 
133c4762a1bSJed Brown   PetscFunctionBegin;
1349566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1359566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Udot, &udot));
136c4762a1bSJed Brown   Pmax    = ctx->Pmax;
1379371c9d4SSatish Balay   J[0][0] = a;
1389371c9d4SSatish Balay   J[0][1] = -ctx->omega_b;
1399371c9d4SSatish Balay   J[1][1] = 2.0 * ctx->H / ctx->omega_s * a + ctx->D;
1409371c9d4SSatish Balay   J[1][0] = Pmax * PetscCosScalar(u[0]);
141c4762a1bSJed Brown 
1429566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES));
1439566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1449566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Udot, &udot));
145c4762a1bSJed Brown 
1469566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1479566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
148c4762a1bSJed Brown   if (A != B) {
1499566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
1509566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
151c4762a1bSJed Brown   }
1523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
153c4762a1bSJed Brown }
154c4762a1bSJed Brown 
155d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSJacobianP(TS ts, PetscReal t, Vec X, Mat A, void *ctx0)
156d71ae5a4SJacob Faibussowitsch {
157c4762a1bSJed Brown   PetscInt     row[] = {0, 1}, col[] = {0};
158c4762a1bSJed Brown   PetscScalar *x, J[2][1];
159c4762a1bSJed Brown   AppCtx      *ctx = (AppCtx *)ctx0;
160c4762a1bSJed Brown 
161c4762a1bSJed Brown   PetscFunctionBeginUser;
1629566063dSJacob Faibussowitsch   PetscCall(VecGetArray(X, &x));
163c4762a1bSJed Brown   J[0][0] = 0;
164c4762a1bSJed Brown   J[1][0] = ctx->omega_s / (2.0 * ctx->H);
1659566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 2, row, 1, col, &J[0][0], INSERT_VALUES));
166c4762a1bSJed Brown 
1679566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1689566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
1693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
170c4762a1bSJed Brown }
171c4762a1bSJed Brown 
172d71ae5a4SJacob Faibussowitsch PetscErrorCode CostIntegrand(TS ts, PetscReal t, Vec U, Vec R, AppCtx *ctx)
173d71ae5a4SJacob Faibussowitsch {
174c4762a1bSJed Brown   PetscScalar       *r;
175c4762a1bSJed Brown   const PetscScalar *u;
176c4762a1bSJed Brown 
177c4762a1bSJed Brown   PetscFunctionBegin;
1789566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1799566063dSJacob Faibussowitsch   PetscCall(VecGetArray(R, &r));
1802f613bf5SBarry Smith   r[0] = ctx->c * PetscPowScalarInt(PetscMax(0., u[0] - ctx->u_s), ctx->beta);
1819566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(R, &r));
1829566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
184c4762a1bSJed Brown }
185c4762a1bSJed Brown 
186c4762a1bSJed Brown /* Transpose of DRDU */
187d71ae5a4SJacob Faibussowitsch PetscErrorCode DRDUJacobianTranspose(TS ts, PetscReal t, Vec U, Mat DRDU, Mat B, AppCtx *ctx)
188d71ae5a4SJacob Faibussowitsch {
189c4762a1bSJed Brown   PetscScalar        ru[2];
190c4762a1bSJed Brown   PetscInt           row[] = {0, 1}, col[] = {0};
191c4762a1bSJed Brown   const PetscScalar *u;
192c4762a1bSJed Brown 
193c4762a1bSJed Brown   PetscFunctionBegin;
1949566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1952f613bf5SBarry Smith   ru[0] = ctx->c * ctx->beta * PetscPowScalarInt(PetscMax(0., u[0] - ctx->u_s), ctx->beta - 1);
196c4762a1bSJed Brown   ru[1] = 0;
1979566063dSJacob Faibussowitsch   PetscCall(MatSetValues(DRDU, 2, row, 1, col, ru, INSERT_VALUES));
1989566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1999566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(DRDU, MAT_FINAL_ASSEMBLY));
2009566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(DRDU, MAT_FINAL_ASSEMBLY));
2013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
202c4762a1bSJed Brown }
203c4762a1bSJed Brown 
204d71ae5a4SJacob Faibussowitsch PetscErrorCode DRDPJacobianTranspose(TS ts, PetscReal t, Vec U, Mat DRDP, void *ctx)
205d71ae5a4SJacob Faibussowitsch {
206c4762a1bSJed Brown   PetscFunctionBegin;
2079566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(DRDP));
2089566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(DRDP, MAT_FINAL_ASSEMBLY));
2099566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(DRDP, MAT_FINAL_ASSEMBLY));
2103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
211c4762a1bSJed Brown }
212c4762a1bSJed Brown 
213d71ae5a4SJacob Faibussowitsch PetscErrorCode ComputeSensiP(Vec lambda, Vec mu, AppCtx *ctx)
214d71ae5a4SJacob Faibussowitsch {
215c4762a1bSJed Brown   PetscScalar       *y, sensip;
216c4762a1bSJed Brown   const PetscScalar *x;
217c4762a1bSJed Brown 
218c4762a1bSJed Brown   PetscFunctionBegin;
2199566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(lambda, &x));
2209566063dSJacob Faibussowitsch   PetscCall(VecGetArray(mu, &y));
221c4762a1bSJed Brown   sensip = 1. / PetscSqrtScalar(1. - (ctx->Pm / ctx->Pmax) * (ctx->Pm / ctx->Pmax)) / ctx->Pmax * x[0] + y[0];
222c4762a1bSJed Brown   y[0]   = sensip;
2239566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(mu, &y));
2249566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(lambda, &x));
2253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
226c4762a1bSJed Brown }
227