xref: /petsc/src/ts/tutorials/power_grid/stability_9bus/ex9busadj.c (revision 58d68138c660dfb4e9f5b03334792cd4f2ffd7cc)
1 
2 static char help[] = "Sensitivity analysis applied in power grid stability analysis of WECC 9 bus system.\n\
3 This example is based on the 9-bus (node) example given in the book Power\n\
4 Systems Dynamics and Stability (Chapter 7) by P. Sauer and M. A. Pai.\n\
5 The power grid in this example consists of 9 buses (nodes), 3 generators,\n\
6 3 loads, and 9 transmission lines. The network equations are written\n\
7 in current balance form using rectangular coordinates.\n\n";
8 
9 /*
10    The equations for the stability analysis are described by the DAE. See ex9bus.c for details.
11    The system has 'jumps' due to faults, thus the time interval is split into multiple sections, and TSSolve is called for each of them. But TSAdjointSolve only needs to be called once since the whole trajectory has been saved in the forward run.
12    The code computes the sensitivity of a final state w.r.t. initial conditions.
13 */
14 
15 #include <petscts.h>
16 #include <petscdm.h>
17 #include <petscdmda.h>
18 #include <petscdmcomposite.h>
19 
20 #define freq 60
21 #define w_s  (2 * PETSC_PI * freq)
22 
23 /* Sizes and indices */
24 const PetscInt nbus    = 9;         /* Number of network buses */
25 const PetscInt ngen    = 3;         /* Number of generators */
26 const PetscInt nload   = 3;         /* Number of loads */
27 const PetscInt gbus[3] = {0, 1, 2}; /* Buses at which generators are incident */
28 const PetscInt lbus[3] = {4, 5, 7}; /* Buses at which loads are incident */
29 
30 /* Generator real and reactive powers (found via loadflow) */
31 const PetscScalar PG[3]   = {0.716786142395021, 1.630000000000000, 0.850000000000000};
32 const PetscScalar QG[3]   = {0.270702180178785, 0.066120127797275, -0.108402221791588};
33 /* Generator constants */
34 const PetscScalar H[3]    = {23.64, 6.4, 3.01};       /* Inertia constant */
35 const PetscScalar Rs[3]   = {0.0, 0.0, 0.0};          /* Stator Resistance */
36 const PetscScalar Xd[3]   = {0.146, 0.8958, 1.3125};  /* d-axis reactance */
37 const PetscScalar Xdp[3]  = {0.0608, 0.1198, 0.1813}; /* d-axis transient reactance */
38 const PetscScalar Xq[3]   = {0.4360, 0.8645, 1.2578}; /* q-axis reactance Xq(1) set to 0.4360, value given in text 0.0969 */
39 const PetscScalar Xqp[3]  = {0.0969, 0.1969, 0.25};   /* q-axis transient reactance */
40 const PetscScalar Td0p[3] = {8.96, 6.0, 5.89};        /* d-axis open circuit time constant */
41 const PetscScalar Tq0p[3] = {0.31, 0.535, 0.6};       /* q-axis open circuit time constant */
42 PetscScalar       M[3];                               /* M = 2*H/w_s */
43 PetscScalar       D[3];                               /* D = 0.1*M */
44 
45 PetscScalar       TM[3]; /* Mechanical Torque */
46 /* Exciter system constants */
47 const PetscScalar KA[3] = {20.0, 20.0, 20.0};    /* Voltage regulartor gain constant */
48 const PetscScalar TA[3] = {0.2, 0.2, 0.2};       /* Voltage regulator time constant */
49 const PetscScalar KE[3] = {1.0, 1.0, 1.0};       /* Exciter gain constant */
50 const PetscScalar TE[3] = {0.314, 0.314, 0.314}; /* Exciter time constant */
51 const PetscScalar KF[3] = {0.063, 0.063, 0.063}; /* Feedback stabilizer gain constant */
52 const PetscScalar TF[3] = {0.35, 0.35, 0.35};    /* Feedback stabilizer time constant */
53 const PetscScalar k1[3] = {0.0039, 0.0039, 0.0039};
54 const PetscScalar k2[3] = {1.555, 1.555, 1.555}; /* k1 and k2 for calculating the saturation function SE = k1*exp(k2*Efd) */
55 
56 PetscScalar       Vref[3];
57 /* Load constants
58   We use a composite load model that describes the load and reactive powers at each time instant as follows
59   P(t) = \sum\limits_{i=0}^ld_nsegsp \ld_alphap_i*P_D0(\frac{V_m(t)}{V_m0})^\ld_betap_i
60   Q(t) = \sum\limits_{i=0}^ld_nsegsq \ld_alphaq_i*Q_D0(\frac{V_m(t)}{V_m0})^\ld_betaq_i
61   where
62     ld_nsegsp,ld_nsegsq - Number of individual load models for real and reactive power loads
63     ld_alphap,ld_alphap - Percentage contribution (weights) or loads
64     P_D0                - Real power load
65     Q_D0                - Reactive power load
66     V_m(t)              - Voltage magnitude at time t
67     V_m0                - Voltage magnitude at t = 0
68     ld_betap, ld_betaq  - exponents describing the load model for real and reactive part
69 
70     Note: All loads have the same characteristic currently.
71 */
72 const PetscScalar PD0[3]       = {1.25, 0.9, 1.0};
73 const PetscScalar QD0[3]       = {0.5, 0.3, 0.35};
74 const PetscInt    ld_nsegsp[3] = {3, 3, 3};
75 const PetscScalar ld_alphap[3] = {1.0, 0.0, 0.0};
76 const PetscScalar ld_betap[3]  = {2.0, 1.0, 0.0};
77 const PetscInt    ld_nsegsq[3] = {3, 3, 3};
78 const PetscScalar ld_alphaq[3] = {1.0, 0.0, 0.0};
79 const PetscScalar ld_betaq[3]  = {2.0, 1.0, 0.0};
80 
81 typedef struct {
82   DM          dmgen, dmnet;        /* DMs to manage generator and network subsystem */
83   DM          dmpgrid;             /* Composite DM to manage the entire power grid */
84   Mat         Ybus;                /* Network admittance matrix */
85   Vec         V0;                  /* Initial voltage vector (Power flow solution) */
86   PetscReal   tfaulton, tfaultoff; /* Fault on and off times */
87   PetscInt    faultbus;            /* Fault bus */
88   PetscScalar Rfault;
89   PetscReal   t0, tmax;
90   PetscInt    neqs_gen, neqs_net, neqs_pgrid;
91   PetscBool   alg_flg;
92   PetscReal   t;
93   IS          is_diff; /* indices for differential equations */
94   IS          is_alg;  /* indices for algebraic equations */
95 } Userctx;
96 
97 /* Converts from machine frame (dq) to network (phase a real,imag) reference frame */
98 PetscErrorCode dq2ri(PetscScalar Fd, PetscScalar Fq, PetscScalar delta, PetscScalar *Fr, PetscScalar *Fi) {
99   PetscFunctionBegin;
100   *Fr = Fd * PetscSinScalar(delta) + Fq * PetscCosScalar(delta);
101   *Fi = -Fd * PetscCosScalar(delta) + Fq * PetscSinScalar(delta);
102   PetscFunctionReturn(0);
103 }
104 
105 /* Converts from network frame ([phase a real,imag) to machine (dq) reference frame */
106 PetscErrorCode ri2dq(PetscScalar Fr, PetscScalar Fi, PetscScalar delta, PetscScalar *Fd, PetscScalar *Fq) {
107   PetscFunctionBegin;
108   *Fd = Fr * PetscSinScalar(delta) - Fi * PetscCosScalar(delta);
109   *Fq = Fr * PetscCosScalar(delta) + Fi * PetscSinScalar(delta);
110   PetscFunctionReturn(0);
111 }
112 
113 PetscErrorCode SetInitialGuess(Vec X, Userctx *user) {
114   Vec          Xgen, Xnet;
115   PetscScalar *xgen, *xnet;
116   PetscInt     i, idx = 0;
117   PetscScalar  Vr, Vi, IGr, IGi, Vm, Vm2;
118   PetscScalar  Eqp, Edp, delta;
119   PetscScalar  Efd, RF, VR; /* Exciter variables */
120   PetscScalar  Id, Iq;      /* Generator dq axis currents */
121   PetscScalar  theta, Vd, Vq, SE;
122 
123   PetscFunctionBegin;
124   M[0] = 2 * H[0] / w_s;
125   M[1] = 2 * H[1] / w_s;
126   M[2] = 2 * H[2] / w_s;
127   D[0] = 0.1 * M[0];
128   D[1] = 0.1 * M[1];
129   D[2] = 0.1 * M[2];
130 
131   PetscCall(DMCompositeGetLocalVectors(user->dmpgrid, &Xgen, &Xnet));
132 
133   /* Network subsystem initialization */
134   PetscCall(VecCopy(user->V0, Xnet));
135 
136   /* Generator subsystem initialization */
137   PetscCall(VecGetArray(Xgen, &xgen));
138   PetscCall(VecGetArray(Xnet, &xnet));
139 
140   for (i = 0; i < ngen; i++) {
141     Vr  = xnet[2 * gbus[i]];     /* Real part of generator terminal voltage */
142     Vi  = xnet[2 * gbus[i] + 1]; /* Imaginary part of the generator terminal voltage */
143     Vm  = PetscSqrtScalar(Vr * Vr + Vi * Vi);
144     Vm2 = Vm * Vm;
145     IGr = (Vr * PG[i] + Vi * QG[i]) / Vm2;
146     IGi = (Vi * PG[i] - Vr * QG[i]) / Vm2;
147 
148     delta = PetscAtan2Real(Vi + Xq[i] * IGr, Vr - Xq[i] * IGi); /* Machine angle */
149 
150     theta = PETSC_PI / 2.0 - delta;
151 
152     Id = IGr * PetscCosScalar(theta) - IGi * PetscSinScalar(theta); /* d-axis stator current */
153     Iq = IGr * PetscSinScalar(theta) + IGi * PetscCosScalar(theta); /* q-axis stator current */
154 
155     Vd = Vr * PetscCosScalar(theta) - Vi * PetscSinScalar(theta);
156     Vq = Vr * PetscSinScalar(theta) + Vi * PetscCosScalar(theta);
157 
158     Edp = Vd + Rs[i] * Id - Xqp[i] * Iq; /* d-axis transient EMF */
159     Eqp = Vq + Rs[i] * Iq + Xdp[i] * Id; /* q-axis transient EMF */
160 
161     TM[i] = PG[i];
162 
163     /* The generator variables are ordered as [Eqp,Edp,delta,w,Id,Iq] */
164     xgen[idx]     = Eqp;
165     xgen[idx + 1] = Edp;
166     xgen[idx + 2] = delta;
167     xgen[idx + 3] = w_s;
168 
169     idx = idx + 4;
170 
171     xgen[idx]     = Id;
172     xgen[idx + 1] = Iq;
173 
174     idx = idx + 2;
175 
176     /* Exciter */
177     Efd = Eqp + (Xd[i] - Xdp[i]) * Id;
178     SE  = k1[i] * PetscExpScalar(k2[i] * Efd);
179     VR  = KE[i] * Efd + SE;
180     RF  = KF[i] * Efd / TF[i];
181 
182     xgen[idx]     = Efd;
183     xgen[idx + 1] = RF;
184     xgen[idx + 2] = VR;
185 
186     Vref[i] = Vm + (VR / KA[i]);
187 
188     idx = idx + 3;
189   }
190 
191   PetscCall(VecRestoreArray(Xgen, &xgen));
192   PetscCall(VecRestoreArray(Xnet, &xnet));
193 
194   /* PetscCall(VecView(Xgen,0)); */
195   PetscCall(DMCompositeGather(user->dmpgrid, INSERT_VALUES, X, Xgen, Xnet));
196   PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid, &Xgen, &Xnet));
197   PetscFunctionReturn(0);
198 }
199 
200 /* Computes F = [f(x,y);g(x,y)] */
201 PetscErrorCode ResidualFunction(SNES snes, Vec X, Vec F, Userctx *user) {
202   Vec          Xgen, Xnet, Fgen, Fnet;
203   PetscScalar *xgen, *xnet, *fgen, *fnet;
204   PetscInt     i, idx = 0;
205   PetscScalar  Vr, Vi, Vm, Vm2;
206   PetscScalar  Eqp, Edp, delta, w; /* Generator variables */
207   PetscScalar  Efd, RF, VR;        /* Exciter variables */
208   PetscScalar  Id, Iq;             /* Generator dq axis currents */
209   PetscScalar  Vd, Vq, SE;
210   PetscScalar  IGr, IGi, IDr, IDi;
211   PetscScalar  Zdq_inv[4], det;
212   PetscScalar  PD, QD, Vm0, *v0;
213   PetscInt     k;
214 
215   PetscFunctionBegin;
216   PetscCall(VecZeroEntries(F));
217   PetscCall(DMCompositeGetLocalVectors(user->dmpgrid, &Xgen, &Xnet));
218   PetscCall(DMCompositeGetLocalVectors(user->dmpgrid, &Fgen, &Fnet));
219   PetscCall(DMCompositeScatter(user->dmpgrid, X, Xgen, Xnet));
220   PetscCall(DMCompositeScatter(user->dmpgrid, F, Fgen, Fnet));
221 
222   /* Network current balance residual IG + Y*V + IL = 0. Only YV is added here.
223      The generator current injection, IG, and load current injection, ID are added later
224   */
225   /* Note that the values in Ybus are stored assuming the imaginary current balance
226      equation is ordered first followed by real current balance equation for each bus.
227      Thus imaginary current contribution goes in location 2*i, and
228      real current contribution in 2*i+1
229   */
230   PetscCall(MatMult(user->Ybus, Xnet, Fnet));
231 
232   PetscCall(VecGetArray(Xgen, &xgen));
233   PetscCall(VecGetArray(Xnet, &xnet));
234   PetscCall(VecGetArray(Fgen, &fgen));
235   PetscCall(VecGetArray(Fnet, &fnet));
236 
237   /* Generator subsystem */
238   for (i = 0; i < ngen; i++) {
239     Eqp   = xgen[idx];
240     Edp   = xgen[idx + 1];
241     delta = xgen[idx + 2];
242     w     = xgen[idx + 3];
243     Id    = xgen[idx + 4];
244     Iq    = xgen[idx + 5];
245     Efd   = xgen[idx + 6];
246     RF    = xgen[idx + 7];
247     VR    = xgen[idx + 8];
248 
249     /* Generator differential equations */
250     fgen[idx]     = (Eqp + (Xd[i] - Xdp[i]) * Id - Efd) / Td0p[i];
251     fgen[idx + 1] = (Edp - (Xq[i] - Xqp[i]) * Iq) / Tq0p[i];
252     fgen[idx + 2] = -w + w_s;
253     fgen[idx + 3] = (-TM[i] + Edp * Id + Eqp * Iq + (Xqp[i] - Xdp[i]) * Id * Iq + D[i] * (w - w_s)) / M[i];
254 
255     Vr = xnet[2 * gbus[i]];     /* Real part of generator terminal voltage */
256     Vi = xnet[2 * gbus[i] + 1]; /* Imaginary part of the generator terminal voltage */
257 
258     PetscCall(ri2dq(Vr, Vi, delta, &Vd, &Vq));
259     /* Algebraic equations for stator currents */
260 
261     det = Rs[i] * Rs[i] + Xdp[i] * Xqp[i];
262 
263     Zdq_inv[0] = Rs[i] / det;
264     Zdq_inv[1] = Xqp[i] / det;
265     Zdq_inv[2] = -Xdp[i] / det;
266     Zdq_inv[3] = Rs[i] / det;
267 
268     fgen[idx + 4] = Zdq_inv[0] * (-Edp + Vd) + Zdq_inv[1] * (-Eqp + Vq) + Id;
269     fgen[idx + 5] = Zdq_inv[2] * (-Edp + Vd) + Zdq_inv[3] * (-Eqp + Vq) + Iq;
270 
271     /* Add generator current injection to network */
272     PetscCall(dq2ri(Id, Iq, delta, &IGr, &IGi));
273 
274     fnet[2 * gbus[i]] -= IGi;
275     fnet[2 * gbus[i] + 1] -= IGr;
276 
277     Vm = PetscSqrtScalar(Vd * Vd + Vq * Vq);
278 
279     SE = k1[i] * PetscExpScalar(k2[i] * Efd);
280 
281     /* Exciter differential equations */
282     fgen[idx + 6] = (KE[i] * Efd + SE - VR) / TE[i];
283     fgen[idx + 7] = (RF - KF[i] * Efd / TF[i]) / TF[i];
284     fgen[idx + 8] = (VR - KA[i] * RF + KA[i] * KF[i] * Efd / TF[i] - KA[i] * (Vref[i] - Vm)) / TA[i];
285 
286     idx = idx + 9;
287   }
288 
289   PetscCall(VecGetArray(user->V0, &v0));
290   for (i = 0; i < nload; i++) {
291     Vr  = xnet[2 * lbus[i]];     /* Real part of load bus voltage */
292     Vi  = xnet[2 * lbus[i] + 1]; /* Imaginary part of the load bus voltage */
293     Vm  = PetscSqrtScalar(Vr * Vr + Vi * Vi);
294     Vm2 = Vm * Vm;
295     Vm0 = PetscSqrtScalar(v0[2 * lbus[i]] * v0[2 * lbus[i]] + v0[2 * lbus[i] + 1] * v0[2 * lbus[i] + 1]);
296     PD = QD = 0.0;
297     for (k = 0; k < ld_nsegsp[i]; k++) PD += ld_alphap[k] * PD0[i] * PetscPowScalar((Vm / Vm0), ld_betap[k]);
298     for (k = 0; k < ld_nsegsq[i]; k++) QD += ld_alphaq[k] * QD0[i] * PetscPowScalar((Vm / Vm0), ld_betaq[k]);
299 
300     /* Load currents */
301     IDr = (PD * Vr + QD * Vi) / Vm2;
302     IDi = (-QD * Vr + PD * Vi) / Vm2;
303 
304     fnet[2 * lbus[i]] += IDi;
305     fnet[2 * lbus[i] + 1] += IDr;
306   }
307   PetscCall(VecRestoreArray(user->V0, &v0));
308 
309   PetscCall(VecRestoreArray(Xgen, &xgen));
310   PetscCall(VecRestoreArray(Xnet, &xnet));
311   PetscCall(VecRestoreArray(Fgen, &fgen));
312   PetscCall(VecRestoreArray(Fnet, &fnet));
313 
314   PetscCall(DMCompositeGather(user->dmpgrid, INSERT_VALUES, F, Fgen, Fnet));
315   PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid, &Xgen, &Xnet));
316   PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid, &Fgen, &Fnet));
317   PetscFunctionReturn(0);
318 }
319 
320 /* \dot{x} - f(x,y)
321      g(x,y) = 0
322  */
323 PetscErrorCode IFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, Userctx *user) {
324   SNES               snes;
325   PetscScalar       *f;
326   const PetscScalar *xdot;
327   PetscInt           i;
328 
329   PetscFunctionBegin;
330   user->t = t;
331 
332   PetscCall(TSGetSNES(ts, &snes));
333   PetscCall(ResidualFunction(snes, X, F, user));
334   PetscCall(VecGetArray(F, &f));
335   PetscCall(VecGetArrayRead(Xdot, &xdot));
336   for (i = 0; i < ngen; i++) {
337     f[9 * i] += xdot[9 * i];
338     f[9 * i + 1] += xdot[9 * i + 1];
339     f[9 * i + 2] += xdot[9 * i + 2];
340     f[9 * i + 3] += xdot[9 * i + 3];
341     f[9 * i + 6] += xdot[9 * i + 6];
342     f[9 * i + 7] += xdot[9 * i + 7];
343     f[9 * i + 8] += xdot[9 * i + 8];
344   }
345   PetscCall(VecRestoreArray(F, &f));
346   PetscCall(VecRestoreArrayRead(Xdot, &xdot));
347   PetscFunctionReturn(0);
348 }
349 
350 /* This function is used for solving the algebraic system only during fault on and
351    off times. It computes the entire F and then zeros out the part corresponding to
352    differential equations
353  F = [0;g(y)];
354 */
355 PetscErrorCode AlgFunction(SNES snes, Vec X, Vec F, void *ctx) {
356   Userctx     *user = (Userctx *)ctx;
357   PetscScalar *f;
358   PetscInt     i;
359 
360   PetscFunctionBegin;
361   PetscCall(ResidualFunction(snes, X, F, user));
362   PetscCall(VecGetArray(F, &f));
363   for (i = 0; i < ngen; i++) {
364     f[9 * i]     = 0;
365     f[9 * i + 1] = 0;
366     f[9 * i + 2] = 0;
367     f[9 * i + 3] = 0;
368     f[9 * i + 6] = 0;
369     f[9 * i + 7] = 0;
370     f[9 * i + 8] = 0;
371   }
372   PetscCall(VecRestoreArray(F, &f));
373   PetscFunctionReturn(0);
374 }
375 
376 PetscErrorCode PreallocateJacobian(Mat J, Userctx *user) {
377   PetscInt *d_nnz;
378   PetscInt  i, idx = 0, start = 0;
379   PetscInt  ncols;
380 
381   PetscFunctionBegin;
382   PetscCall(PetscMalloc1(user->neqs_pgrid, &d_nnz));
383   for (i = 0; i < user->neqs_pgrid; i++) d_nnz[i] = 0;
384   /* Generator subsystem */
385   for (i = 0; i < ngen; i++) {
386     d_nnz[idx] += 3;
387     d_nnz[idx + 1] += 2;
388     d_nnz[idx + 2] += 2;
389     d_nnz[idx + 3] += 5;
390     d_nnz[idx + 4] += 6;
391     d_nnz[idx + 5] += 6;
392 
393     d_nnz[user->neqs_gen + 2 * gbus[i]] += 3;
394     d_nnz[user->neqs_gen + 2 * gbus[i] + 1] += 3;
395 
396     d_nnz[idx + 6] += 2;
397     d_nnz[idx + 7] += 2;
398     d_nnz[idx + 8] += 5;
399 
400     idx = idx + 9;
401   }
402 
403   start = user->neqs_gen;
404 
405   for (i = 0; i < nbus; i++) {
406     PetscCall(MatGetRow(user->Ybus, 2 * i, &ncols, NULL, NULL));
407     d_nnz[start + 2 * i] += ncols;
408     d_nnz[start + 2 * i + 1] += ncols;
409     PetscCall(MatRestoreRow(user->Ybus, 2 * i, &ncols, NULL, NULL));
410   }
411 
412   PetscCall(MatSeqAIJSetPreallocation(J, 0, d_nnz));
413 
414   PetscCall(PetscFree(d_nnz));
415   PetscFunctionReturn(0);
416 }
417 
418 /*
419    J = [-df_dx, -df_dy
420         dg_dx, dg_dy]
421 */
422 PetscErrorCode ResidualJacobian(SNES snes, Vec X, Mat J, Mat B, void *ctx) {
423   Userctx           *user = (Userctx *)ctx;
424   Vec                Xgen, Xnet;
425   PetscScalar       *xgen, *xnet;
426   PetscInt           i, idx = 0;
427   PetscScalar        Vr, Vi, Vm, Vm2;
428   PetscScalar        Eqp, Edp, delta; /* Generator variables */
429   PetscScalar        Efd;             /* Exciter variables */
430   PetscScalar        Id, Iq;          /* Generator dq axis currents */
431   PetscScalar        Vd, Vq;
432   PetscScalar        val[10];
433   PetscInt           row[2], col[10];
434   PetscInt           net_start = user->neqs_gen;
435   PetscScalar        Zdq_inv[4], det;
436   PetscScalar        dVd_dVr, dVd_dVi, dVq_dVr, dVq_dVi, dVd_ddelta, dVq_ddelta;
437   PetscScalar        dIGr_ddelta, dIGi_ddelta, dIGr_dId, dIGr_dIq, dIGi_dId, dIGi_dIq;
438   PetscScalar        dSE_dEfd;
439   PetscScalar        dVm_dVd, dVm_dVq, dVm_dVr, dVm_dVi;
440   PetscInt           ncols;
441   const PetscInt    *cols;
442   const PetscScalar *yvals;
443   PetscInt           k;
444   PetscScalar        PD, QD, Vm0, *v0, Vm4;
445   PetscScalar        dPD_dVr, dPD_dVi, dQD_dVr, dQD_dVi;
446   PetscScalar        dIDr_dVr, dIDr_dVi, dIDi_dVr, dIDi_dVi;
447 
448   PetscFunctionBegin;
449   PetscCall(MatZeroEntries(B));
450   PetscCall(DMCompositeGetLocalVectors(user->dmpgrid, &Xgen, &Xnet));
451   PetscCall(DMCompositeScatter(user->dmpgrid, X, Xgen, Xnet));
452 
453   PetscCall(VecGetArray(Xgen, &xgen));
454   PetscCall(VecGetArray(Xnet, &xnet));
455 
456   /* Generator subsystem */
457   for (i = 0; i < ngen; i++) {
458     Eqp   = xgen[idx];
459     Edp   = xgen[idx + 1];
460     delta = xgen[idx + 2];
461     Id    = xgen[idx + 4];
462     Iq    = xgen[idx + 5];
463     Efd   = xgen[idx + 6];
464 
465     /*    fgen[idx]   = (Eqp + (Xd[i] - Xdp[i])*Id - Efd)/Td0p[i]; */
466     row[0] = idx;
467     col[0] = idx;
468     col[1] = idx + 4;
469     col[2] = idx + 6;
470     val[0] = 1 / Td0p[i];
471     val[1] = (Xd[i] - Xdp[i]) / Td0p[i];
472     val[2] = -1 / Td0p[i];
473 
474     PetscCall(MatSetValues(J, 1, row, 3, col, val, INSERT_VALUES));
475 
476     /*    fgen[idx+1] = (Edp - (Xq[i] - Xqp[i])*Iq)/Tq0p[i]; */
477     row[0] = idx + 1;
478     col[0] = idx + 1;
479     col[1] = idx + 5;
480     val[0] = 1 / Tq0p[i];
481     val[1] = -(Xq[i] - Xqp[i]) / Tq0p[i];
482     PetscCall(MatSetValues(J, 1, row, 2, col, val, INSERT_VALUES));
483 
484     /*    fgen[idx+2] = - w + w_s; */
485     row[0] = idx + 2;
486     col[0] = idx + 2;
487     col[1] = idx + 3;
488     val[0] = 0;
489     val[1] = -1;
490     PetscCall(MatSetValues(J, 1, row, 2, col, val, INSERT_VALUES));
491 
492     /*    fgen[idx+3] = (-TM[i] + Edp*Id + Eqp*Iq + (Xqp[i] - Xdp[i])*Id*Iq + D[i]*(w - w_s))/M[i]; */
493     row[0] = idx + 3;
494     col[0] = idx;
495     col[1] = idx + 1;
496     col[2] = idx + 3;
497     col[3] = idx + 4;
498     col[4] = idx + 5;
499     val[0] = Iq / M[i];
500     val[1] = Id / M[i];
501     val[2] = D[i] / M[i];
502     val[3] = (Edp + (Xqp[i] - Xdp[i]) * Iq) / M[i];
503     val[4] = (Eqp + (Xqp[i] - Xdp[i]) * Id) / M[i];
504     PetscCall(MatSetValues(J, 1, row, 5, col, val, INSERT_VALUES));
505 
506     Vr = xnet[2 * gbus[i]];     /* Real part of generator terminal voltage */
507     Vi = xnet[2 * gbus[i] + 1]; /* Imaginary part of the generator terminal voltage */
508     PetscCall(ri2dq(Vr, Vi, delta, &Vd, &Vq));
509 
510     det = Rs[i] * Rs[i] + Xdp[i] * Xqp[i];
511 
512     Zdq_inv[0] = Rs[i] / det;
513     Zdq_inv[1] = Xqp[i] / det;
514     Zdq_inv[2] = -Xdp[i] / det;
515     Zdq_inv[3] = Rs[i] / det;
516 
517     dVd_dVr    = PetscSinScalar(delta);
518     dVd_dVi    = -PetscCosScalar(delta);
519     dVq_dVr    = PetscCosScalar(delta);
520     dVq_dVi    = PetscSinScalar(delta);
521     dVd_ddelta = Vr * PetscCosScalar(delta) + Vi * PetscSinScalar(delta);
522     dVq_ddelta = -Vr * PetscSinScalar(delta) + Vi * PetscCosScalar(delta);
523 
524     /*    fgen[idx+4] = Zdq_inv[0]*(-Edp + Vd) + Zdq_inv[1]*(-Eqp + Vq) + Id; */
525     row[0] = idx + 4;
526     col[0] = idx;
527     col[1] = idx + 1;
528     col[2] = idx + 2;
529     val[0] = -Zdq_inv[1];
530     val[1] = -Zdq_inv[0];
531     val[2] = Zdq_inv[0] * dVd_ddelta + Zdq_inv[1] * dVq_ddelta;
532     col[3] = idx + 4;
533     col[4] = net_start + 2 * gbus[i];
534     col[5] = net_start + 2 * gbus[i] + 1;
535     val[3] = 1;
536     val[4] = Zdq_inv[0] * dVd_dVr + Zdq_inv[1] * dVq_dVr;
537     val[5] = Zdq_inv[0] * dVd_dVi + Zdq_inv[1] * dVq_dVi;
538     PetscCall(MatSetValues(J, 1, row, 6, col, val, INSERT_VALUES));
539 
540     /*  fgen[idx+5] = Zdq_inv[2]*(-Edp + Vd) + Zdq_inv[3]*(-Eqp + Vq) + Iq; */
541     row[0] = idx + 5;
542     col[0] = idx;
543     col[1] = idx + 1;
544     col[2] = idx + 2;
545     val[0] = -Zdq_inv[3];
546     val[1] = -Zdq_inv[2];
547     val[2] = Zdq_inv[2] * dVd_ddelta + Zdq_inv[3] * dVq_ddelta;
548     col[3] = idx + 5;
549     col[4] = net_start + 2 * gbus[i];
550     col[5] = net_start + 2 * gbus[i] + 1;
551     val[3] = 1;
552     val[4] = Zdq_inv[2] * dVd_dVr + Zdq_inv[3] * dVq_dVr;
553     val[5] = Zdq_inv[2] * dVd_dVi + Zdq_inv[3] * dVq_dVi;
554     PetscCall(MatSetValues(J, 1, row, 6, col, val, INSERT_VALUES));
555 
556     dIGr_ddelta = Id * PetscCosScalar(delta) - Iq * PetscSinScalar(delta);
557     dIGi_ddelta = Id * PetscSinScalar(delta) + Iq * PetscCosScalar(delta);
558     dIGr_dId    = PetscSinScalar(delta);
559     dIGr_dIq    = PetscCosScalar(delta);
560     dIGi_dId    = -PetscCosScalar(delta);
561     dIGi_dIq    = PetscSinScalar(delta);
562 
563     /* fnet[2*gbus[i]]   -= IGi; */
564     row[0] = net_start + 2 * gbus[i];
565     col[0] = idx + 2;
566     col[1] = idx + 4;
567     col[2] = idx + 5;
568     val[0] = -dIGi_ddelta;
569     val[1] = -dIGi_dId;
570     val[2] = -dIGi_dIq;
571     PetscCall(MatSetValues(J, 1, row, 3, col, val, INSERT_VALUES));
572 
573     /* fnet[2*gbus[i]+1]   -= IGr; */
574     row[0] = net_start + 2 * gbus[i] + 1;
575     col[0] = idx + 2;
576     col[1] = idx + 4;
577     col[2] = idx + 5;
578     val[0] = -dIGr_ddelta;
579     val[1] = -dIGr_dId;
580     val[2] = -dIGr_dIq;
581     PetscCall(MatSetValues(J, 1, row, 3, col, val, INSERT_VALUES));
582 
583     Vm = PetscSqrtScalar(Vd * Vd + Vq * Vq);
584 
585     /*    fgen[idx+6] = (KE[i]*Efd + SE - VR)/TE[i]; */
586     /*    SE  = k1[i]*PetscExpScalar(k2[i]*Efd); */
587     dSE_dEfd = k1[i] * k2[i] * PetscExpScalar(k2[i] * Efd);
588 
589     row[0] = idx + 6;
590     col[0] = idx + 6;
591     col[1] = idx + 8;
592     val[0] = (KE[i] + dSE_dEfd) / TE[i];
593     val[1] = -1 / TE[i];
594     PetscCall(MatSetValues(J, 1, row, 2, col, val, INSERT_VALUES));
595 
596     /* Exciter differential equations */
597 
598     /*    fgen[idx+7] = (RF - KF[i]*Efd/TF[i])/TF[i]; */
599     row[0] = idx + 7;
600     col[0] = idx + 6;
601     col[1] = idx + 7;
602     val[0] = (-KF[i] / TF[i]) / TF[i];
603     val[1] = 1 / TF[i];
604     PetscCall(MatSetValues(J, 1, row, 2, col, val, INSERT_VALUES));
605 
606     /*    fgen[idx+8] = (VR - KA[i]*RF + KA[i]*KF[i]*Efd/TF[i] - KA[i]*(Vref[i] - Vm))/TA[i]; */
607     /* Vm = (Vd^2 + Vq^2)^0.5; */
608     dVm_dVd = Vd / Vm;
609     dVm_dVq = Vq / Vm;
610     dVm_dVr = dVm_dVd * dVd_dVr + dVm_dVq * dVq_dVr;
611     dVm_dVi = dVm_dVd * dVd_dVi + dVm_dVq * dVq_dVi;
612     row[0]  = idx + 8;
613     col[0]  = idx + 6;
614     col[1]  = idx + 7;
615     col[2]  = idx + 8;
616     val[0]  = (KA[i] * KF[i] / TF[i]) / TA[i];
617     val[1]  = -KA[i] / TA[i];
618     val[2]  = 1 / TA[i];
619     col[3]  = net_start + 2 * gbus[i];
620     col[4]  = net_start + 2 * gbus[i] + 1;
621     val[3]  = KA[i] * dVm_dVr / TA[i];
622     val[4]  = KA[i] * dVm_dVi / TA[i];
623     PetscCall(MatSetValues(J, 1, row, 5, col, val, INSERT_VALUES));
624     idx = idx + 9;
625   }
626 
627   for (i = 0; i < nbus; i++) {
628     PetscCall(MatGetRow(user->Ybus, 2 * i, &ncols, &cols, &yvals));
629     row[0] = net_start + 2 * i;
630     for (k = 0; k < ncols; k++) {
631       col[k] = net_start + cols[k];
632       val[k] = yvals[k];
633     }
634     PetscCall(MatSetValues(J, 1, row, ncols, col, val, INSERT_VALUES));
635     PetscCall(MatRestoreRow(user->Ybus, 2 * i, &ncols, &cols, &yvals));
636 
637     PetscCall(MatGetRow(user->Ybus, 2 * i + 1, &ncols, &cols, &yvals));
638     row[0] = net_start + 2 * i + 1;
639     for (k = 0; k < ncols; k++) {
640       col[k] = net_start + cols[k];
641       val[k] = yvals[k];
642     }
643     PetscCall(MatSetValues(J, 1, row, ncols, col, val, INSERT_VALUES));
644     PetscCall(MatRestoreRow(user->Ybus, 2 * i + 1, &ncols, &cols, &yvals));
645   }
646 
647   PetscCall(MatAssemblyBegin(J, MAT_FLUSH_ASSEMBLY));
648   PetscCall(MatAssemblyEnd(J, MAT_FLUSH_ASSEMBLY));
649 
650   PetscCall(VecGetArray(user->V0, &v0));
651   for (i = 0; i < nload; i++) {
652     Vr  = xnet[2 * lbus[i]];     /* Real part of load bus voltage */
653     Vi  = xnet[2 * lbus[i] + 1]; /* Imaginary part of the load bus voltage */
654     Vm  = PetscSqrtScalar(Vr * Vr + Vi * Vi);
655     Vm2 = Vm * Vm;
656     Vm4 = Vm2 * Vm2;
657     Vm0 = PetscSqrtScalar(v0[2 * lbus[i]] * v0[2 * lbus[i]] + v0[2 * lbus[i] + 1] * v0[2 * lbus[i] + 1]);
658     PD = QD = 0.0;
659     dPD_dVr = dPD_dVi = dQD_dVr = dQD_dVi = 0.0;
660     for (k = 0; k < ld_nsegsp[i]; k++) {
661       PD += ld_alphap[k] * PD0[i] * PetscPowScalar((Vm / Vm0), ld_betap[k]);
662       dPD_dVr += ld_alphap[k] * ld_betap[k] * PD0[i] * PetscPowScalar((1 / Vm0), ld_betap[k]) * Vr * PetscPowScalar(Vm, (ld_betap[k] - 2));
663       dPD_dVi += ld_alphap[k] * ld_betap[k] * PD0[i] * PetscPowScalar((1 / Vm0), ld_betap[k]) * Vi * PetscPowScalar(Vm, (ld_betap[k] - 2));
664     }
665     for (k = 0; k < ld_nsegsq[i]; k++) {
666       QD += ld_alphaq[k] * QD0[i] * PetscPowScalar((Vm / Vm0), ld_betaq[k]);
667       dQD_dVr += ld_alphaq[k] * ld_betaq[k] * QD0[i] * PetscPowScalar((1 / Vm0), ld_betaq[k]) * Vr * PetscPowScalar(Vm, (ld_betaq[k] - 2));
668       dQD_dVi += ld_alphaq[k] * ld_betaq[k] * QD0[i] * PetscPowScalar((1 / Vm0), ld_betaq[k]) * Vi * PetscPowScalar(Vm, (ld_betaq[k] - 2));
669     }
670 
671     /*    IDr = (PD*Vr + QD*Vi)/Vm2; */
672     /*    IDi = (-QD*Vr + PD*Vi)/Vm2; */
673 
674     dIDr_dVr = (dPD_dVr * Vr + dQD_dVr * Vi + PD) / Vm2 - ((PD * Vr + QD * Vi) * 2 * Vr) / Vm4;
675     dIDr_dVi = (dPD_dVi * Vr + dQD_dVi * Vi + QD) / Vm2 - ((PD * Vr + QD * Vi) * 2 * Vi) / Vm4;
676 
677     dIDi_dVr = (-dQD_dVr * Vr + dPD_dVr * Vi - QD) / Vm2 - ((-QD * Vr + PD * Vi) * 2 * Vr) / Vm4;
678     dIDi_dVi = (-dQD_dVi * Vr + dPD_dVi * Vi + PD) / Vm2 - ((-QD * Vr + PD * Vi) * 2 * Vi) / Vm4;
679 
680     /*    fnet[2*lbus[i]]   += IDi; */
681     row[0] = net_start + 2 * lbus[i];
682     col[0] = net_start + 2 * lbus[i];
683     col[1] = net_start + 2 * lbus[i] + 1;
684     val[0] = dIDi_dVr;
685     val[1] = dIDi_dVi;
686     PetscCall(MatSetValues(J, 1, row, 2, col, val, ADD_VALUES));
687     /*    fnet[2*lbus[i]+1] += IDr; */
688     row[0] = net_start + 2 * lbus[i] + 1;
689     col[0] = net_start + 2 * lbus[i];
690     col[1] = net_start + 2 * lbus[i] + 1;
691     val[0] = dIDr_dVr;
692     val[1] = dIDr_dVi;
693     PetscCall(MatSetValues(J, 1, row, 2, col, val, ADD_VALUES));
694   }
695   PetscCall(VecRestoreArray(user->V0, &v0));
696 
697   PetscCall(VecRestoreArray(Xgen, &xgen));
698   PetscCall(VecRestoreArray(Xnet, &xnet));
699 
700   PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid, &Xgen, &Xnet));
701 
702   PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
703   PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
704   PetscFunctionReturn(0);
705 }
706 
707 /*
708    J = [I, 0
709         dg_dx, dg_dy]
710 */
711 PetscErrorCode AlgJacobian(SNES snes, Vec X, Mat A, Mat B, void *ctx) {
712   Userctx *user = (Userctx *)ctx;
713 
714   PetscFunctionBegin;
715   PetscCall(ResidualJacobian(snes, X, A, B, ctx));
716   PetscCall(MatSetOption(A, MAT_KEEP_NONZERO_PATTERN, PETSC_TRUE));
717   PetscCall(MatZeroRowsIS(A, user->is_diff, 1.0, NULL, NULL));
718   PetscFunctionReturn(0);
719 }
720 
721 /*
722    J = [a*I-df_dx, -df_dy
723         dg_dx, dg_dy]
724 */
725 
726 PetscErrorCode IJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat A, Mat B, Userctx *user) {
727   SNES        snes;
728   PetscScalar atmp = (PetscScalar)a;
729   PetscInt    i, row;
730 
731   PetscFunctionBegin;
732   user->t = t;
733 
734   PetscCall(TSGetSNES(ts, &snes));
735   PetscCall(ResidualJacobian(snes, X, A, B, user));
736   for (i = 0; i < ngen; i++) {
737     row = 9 * i;
738     PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES));
739     row = 9 * i + 1;
740     PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES));
741     row = 9 * i + 2;
742     PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES));
743     row = 9 * i + 3;
744     PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES));
745     row = 9 * i + 6;
746     PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES));
747     row = 9 * i + 7;
748     PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES));
749     row = 9 * i + 8;
750     PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES));
751   }
752   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
753   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
754   PetscFunctionReturn(0);
755 }
756 
757 int main(int argc, char **argv) {
758   TS           ts;
759   SNES         snes_alg;
760   PetscMPIInt  size;
761   Userctx      user;
762   PetscViewer  Xview, Ybusview;
763   Vec          X;
764   Mat          J;
765   PetscInt     i;
766   /* sensitivity context */
767   PetscScalar *y_ptr;
768   Vec          lambda[1];
769   PetscInt    *idx2;
770   Vec          Xdot;
771   Vec          F_alg;
772   PetscInt     row_loc, col_loc;
773   PetscScalar  val;
774 
775   PetscFunctionBeginUser;
776   PetscCall(PetscInitialize(&argc, &argv, "petscoptions", help));
777   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
778   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs");
779 
780   user.neqs_gen   = 9 * ngen; /* # eqs. for generator subsystem */
781   user.neqs_net   = 2 * nbus; /* # eqs. for network subsystem   */
782   user.neqs_pgrid = user.neqs_gen + user.neqs_net;
783 
784   /* Create indices for differential and algebraic equations */
785   PetscCall(PetscMalloc1(7 * ngen, &idx2));
786   for (i = 0; i < ngen; i++) {
787     idx2[7 * i]     = 9 * i;
788     idx2[7 * i + 1] = 9 * i + 1;
789     idx2[7 * i + 2] = 9 * i + 2;
790     idx2[7 * i + 3] = 9 * i + 3;
791     idx2[7 * i + 4] = 9 * i + 6;
792     idx2[7 * i + 5] = 9 * i + 7;
793     idx2[7 * i + 6] = 9 * i + 8;
794   }
795   PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, 7 * ngen, idx2, PETSC_COPY_VALUES, &user.is_diff));
796   PetscCall(ISComplement(user.is_diff, 0, user.neqs_pgrid, &user.is_alg));
797   PetscCall(PetscFree(idx2));
798 
799   /* Read initial voltage vector and Ybus */
800   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "X.bin", FILE_MODE_READ, &Xview));
801   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "Ybus.bin", FILE_MODE_READ, &Ybusview));
802 
803   PetscCall(VecCreate(PETSC_COMM_WORLD, &user.V0));
804   PetscCall(VecSetSizes(user.V0, PETSC_DECIDE, user.neqs_net));
805   PetscCall(VecLoad(user.V0, Xview));
806 
807   PetscCall(MatCreate(PETSC_COMM_WORLD, &user.Ybus));
808   PetscCall(MatSetSizes(user.Ybus, PETSC_DECIDE, PETSC_DECIDE, user.neqs_net, user.neqs_net));
809   PetscCall(MatSetType(user.Ybus, MATBAIJ));
810   /*  PetscCall(MatSetBlockSize(user.Ybus,2)); */
811   PetscCall(MatLoad(user.Ybus, Ybusview));
812 
813   /* Set run time options */
814   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Transient stability fault options", "");
815   {
816     user.tfaulton  = 1.0;
817     user.tfaultoff = 1.2;
818     user.Rfault    = 0.0001;
819     user.faultbus  = 8;
820     PetscCall(PetscOptionsReal("-tfaulton", "", "", user.tfaulton, &user.tfaulton, NULL));
821     PetscCall(PetscOptionsReal("-tfaultoff", "", "", user.tfaultoff, &user.tfaultoff, NULL));
822     PetscCall(PetscOptionsInt("-faultbus", "", "", user.faultbus, &user.faultbus, NULL));
823     user.t0   = 0.0;
824     user.tmax = 5.0;
825     PetscCall(PetscOptionsReal("-t0", "", "", user.t0, &user.t0, NULL));
826     PetscCall(PetscOptionsReal("-tmax", "", "", user.tmax, &user.tmax, NULL));
827   }
828   PetscOptionsEnd();
829 
830   PetscCall(PetscViewerDestroy(&Xview));
831   PetscCall(PetscViewerDestroy(&Ybusview));
832 
833   /* Create DMs for generator and network subsystems */
834   PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, user.neqs_gen, 1, 1, NULL, &user.dmgen));
835   PetscCall(DMSetOptionsPrefix(user.dmgen, "dmgen_"));
836   PetscCall(DMSetFromOptions(user.dmgen));
837   PetscCall(DMSetUp(user.dmgen));
838   PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, user.neqs_net, 1, 1, NULL, &user.dmnet));
839   PetscCall(DMSetOptionsPrefix(user.dmnet, "dmnet_"));
840   PetscCall(DMSetFromOptions(user.dmnet));
841   PetscCall(DMSetUp(user.dmnet));
842   /* Create a composite DM packer and add the two DMs */
843   PetscCall(DMCompositeCreate(PETSC_COMM_WORLD, &user.dmpgrid));
844   PetscCall(DMSetOptionsPrefix(user.dmpgrid, "pgrid_"));
845   PetscCall(DMCompositeAddDM(user.dmpgrid, user.dmgen));
846   PetscCall(DMCompositeAddDM(user.dmpgrid, user.dmnet));
847 
848   PetscCall(DMCreateGlobalVector(user.dmpgrid, &X));
849 
850   PetscCall(MatCreate(PETSC_COMM_WORLD, &J));
851   PetscCall(MatSetSizes(J, PETSC_DECIDE, PETSC_DECIDE, user.neqs_pgrid, user.neqs_pgrid));
852   PetscCall(MatSetFromOptions(J));
853   PetscCall(PreallocateJacobian(J, &user));
854 
855   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
856      Create timestepping solver context
857      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
858   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
859   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
860   PetscCall(TSSetType(ts, TSCN));
861   PetscCall(TSSetIFunction(ts, NULL, (TSIFunction)IFunction, &user));
862   PetscCall(TSSetIJacobian(ts, J, J, (TSIJacobian)IJacobian, &user));
863   PetscCall(TSSetApplicationContext(ts, &user));
864 
865   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
866      Set initial conditions
867    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
868   PetscCall(SetInitialGuess(X, &user));
869   /* Just to set up the Jacobian structure */
870   PetscCall(VecDuplicate(X, &Xdot));
871   PetscCall(IJacobian(ts, 0.0, X, Xdot, 0.0, J, J, &user));
872   PetscCall(VecDestroy(&Xdot));
873 
874   /*
875     Save trajectory of solution so that TSAdjointSolve() may be used
876   */
877   PetscCall(TSSetSaveTrajectory(ts));
878 
879   PetscCall(TSSetMaxTime(ts, user.tmax));
880   PetscCall(TSSetTimeStep(ts, 0.01));
881   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
882   PetscCall(TSSetFromOptions(ts));
883 
884   user.alg_flg = PETSC_FALSE;
885   /* Prefault period */
886   PetscCall(TSSolve(ts, X));
887 
888   /* Create the nonlinear solver for solving the algebraic system */
889   /* Note that although the algebraic system needs to be solved only for
890      Idq and V, we reuse the entire system including xgen. The xgen
891      variables are held constant by setting their residuals to 0 and
892      putting a 1 on the Jacobian diagonal for xgen rows
893   */
894   PetscCall(VecDuplicate(X, &F_alg));
895   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes_alg));
896   PetscCall(SNESSetFunction(snes_alg, F_alg, AlgFunction, &user));
897   PetscCall(MatZeroEntries(J));
898   PetscCall(SNESSetJacobian(snes_alg, J, J, AlgJacobian, &user));
899   PetscCall(SNESSetOptionsPrefix(snes_alg, "alg_"));
900   PetscCall(SNESSetFromOptions(snes_alg));
901 
902   /* Apply disturbance - resistive fault at user.faultbus */
903   /* This is done by adding shunt conductance to the diagonal location
904      in the Ybus matrix */
905   row_loc = 2 * user.faultbus;
906   col_loc = 2 * user.faultbus + 1; /* Location for G */
907   val     = 1 / user.Rfault;
908   PetscCall(MatSetValues(user.Ybus, 1, &row_loc, 1, &col_loc, &val, ADD_VALUES));
909   row_loc = 2 * user.faultbus + 1;
910   col_loc = 2 * user.faultbus; /* Location for G */
911   val     = 1 / user.Rfault;
912   PetscCall(MatSetValues(user.Ybus, 1, &row_loc, 1, &col_loc, &val, ADD_VALUES));
913 
914   PetscCall(MatAssemblyBegin(user.Ybus, MAT_FINAL_ASSEMBLY));
915   PetscCall(MatAssemblyEnd(user.Ybus, MAT_FINAL_ASSEMBLY));
916 
917   user.alg_flg = PETSC_TRUE;
918   /* Solve the algebraic equations */
919   PetscCall(SNESSolve(snes_alg, NULL, X));
920 
921   /* Disturbance period */
922   user.alg_flg = PETSC_FALSE;
923   PetscCall(TSSetTime(ts, user.tfaulton));
924   PetscCall(TSSetMaxTime(ts, user.tfaultoff));
925   PetscCall(TSSolve(ts, X));
926 
927   /* Remove the fault */
928   row_loc = 2 * user.faultbus;
929   col_loc = 2 * user.faultbus + 1;
930   val     = -1 / user.Rfault;
931   PetscCall(MatSetValues(user.Ybus, 1, &row_loc, 1, &col_loc, &val, ADD_VALUES));
932   row_loc = 2 * user.faultbus + 1;
933   col_loc = 2 * user.faultbus;
934   val     = -1 / user.Rfault;
935   PetscCall(MatSetValues(user.Ybus, 1, &row_loc, 1, &col_loc, &val, ADD_VALUES));
936 
937   PetscCall(MatAssemblyBegin(user.Ybus, MAT_FINAL_ASSEMBLY));
938   PetscCall(MatAssemblyEnd(user.Ybus, MAT_FINAL_ASSEMBLY));
939 
940   PetscCall(MatZeroEntries(J));
941 
942   user.alg_flg = PETSC_TRUE;
943 
944   /* Solve the algebraic equations */
945   PetscCall(SNESSolve(snes_alg, NULL, X));
946 
947   /* Post-disturbance period */
948   user.alg_flg = PETSC_TRUE;
949   PetscCall(TSSetTime(ts, user.tfaultoff));
950   PetscCall(TSSetMaxTime(ts, user.tmax));
951   PetscCall(TSSolve(ts, X));
952 
953   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
954      Adjoint model starts here
955      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
956   PetscCall(TSSetPostStep(ts, NULL));
957   PetscCall(MatCreateVecs(J, &lambda[0], NULL));
958   /*   Set initial conditions for the adjoint integration */
959   PetscCall(VecZeroEntries(lambda[0]));
960   PetscCall(VecGetArray(lambda[0], &y_ptr));
961   y_ptr[0] = 1.0;
962   PetscCall(VecRestoreArray(lambda[0], &y_ptr));
963   PetscCall(TSSetCostGradients(ts, 1, lambda, NULL));
964 
965   PetscCall(TSAdjointSolve(ts));
966 
967   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n sensitivity wrt initial conditions: \n"));
968   PetscCall(VecView(lambda[0], PETSC_VIEWER_STDOUT_WORLD));
969   PetscCall(VecDestroy(&lambda[0]));
970 
971   PetscCall(SNESDestroy(&snes_alg));
972   PetscCall(VecDestroy(&F_alg));
973   PetscCall(MatDestroy(&J));
974   PetscCall(MatDestroy(&user.Ybus));
975   PetscCall(VecDestroy(&X));
976   PetscCall(VecDestroy(&user.V0));
977   PetscCall(DMDestroy(&user.dmgen));
978   PetscCall(DMDestroy(&user.dmnet));
979   PetscCall(DMDestroy(&user.dmpgrid));
980   PetscCall(ISDestroy(&user.is_diff));
981   PetscCall(ISDestroy(&user.is_alg));
982   PetscCall(TSDestroy(&ts));
983   PetscCall(PetscFinalize());
984   return 0;
985 }
986 
987 /*TEST
988 
989    build:
990       requires: double !complex !defined(PETSC_USE_64BIT_INDICES)
991 
992    test:
993       args: -viewer_binary_skip_info
994       localrunfiles: petscoptions X.bin Ybus.bin
995 
996    test:
997       suffix: 2
998       args: -viewer_binary_skip_info -ts_type beuler
999       localrunfiles: petscoptions X.bin Ybus.bin
1000 
1001 TEST*/
1002