xref: /petsc/src/ts/tutorials/power_grid/stability_9bus/ex9busopt.c (revision 58d68138c660dfb4e9f5b03334792cd4f2ffd7cc)
1 static char help[] = "Application of adjoint sensitivity analysis for power grid stability analysis of WECC 9 bus system.\n\
2 This example is based on the 9-bus (node) example given in the book Power\n\
3 Systems Dynamics and Stability (Chapter 7) by P. Sauer and M. A. Pai.\n\
4 The power grid in this example consists of 9 buses (nodes), 3 generators,\n\
5 3 loads, and 9 transmission lines. The network equations are written\n\
6 in current balance form using rectangular coordinates.\n\n";
7 
8 /*
9   This code demonstrates how to solve a DAE-constrained optimization problem with TAO, TSAdjoint and TS.
10   The objectivie is to find optimal parameter PG for each generator to minizie the frequency violations due to faults.
11   The problem features discontinuities and a cost function in integral form.
12   The gradient is computed with the discrete adjoint of an implicit theta method, see ex9busadj.c for details.
13 */
14 
15 #include <petsctao.h>
16 #include <petscts.h>
17 #include <petscdm.h>
18 #include <petscdmda.h>
19 #include <petscdmcomposite.h>
20 #include <petsctime.h>
21 
22 PetscErrorCode FormFunctionGradient(Tao, Vec, PetscReal *, Vec, void *);
23 
24 #define freq 60
25 #define w_s  (2 * PETSC_PI * freq)
26 
27 /* Sizes and indices */
28 const PetscInt nbus    = 9;         /* Number of network buses */
29 const PetscInt ngen    = 3;         /* Number of generators */
30 const PetscInt nload   = 3;         /* Number of loads */
31 const PetscInt gbus[3] = {0, 1, 2}; /* Buses at which generators are incident */
32 const PetscInt lbus[3] = {4, 5, 7}; /* Buses at which loads are incident */
33 
34 /* Generator real and reactive powers (found via loadflow) */
35 PetscScalar PG[3] = {0.69, 1.59, 0.69};
36 /* PetscScalar PG[3] = {0.716786142395021,1.630000000000000,0.850000000000000};*/
37 
38 const PetscScalar QG[3]   = {0.270702180178785, 0.066120127797275, -0.108402221791588};
39 /* Generator constants */
40 const PetscScalar H[3]    = {23.64, 6.4, 3.01};       /* Inertia constant */
41 const PetscScalar Rs[3]   = {0.0, 0.0, 0.0};          /* Stator Resistance */
42 const PetscScalar Xd[3]   = {0.146, 0.8958, 1.3125};  /* d-axis reactance */
43 const PetscScalar Xdp[3]  = {0.0608, 0.1198, 0.1813}; /* d-axis transient reactance */
44 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 */
45 const PetscScalar Xqp[3]  = {0.0969, 0.1969, 0.25};   /* q-axis transient reactance */
46 const PetscScalar Td0p[3] = {8.96, 6.0, 5.89};        /* d-axis open circuit time constant */
47 const PetscScalar Tq0p[3] = {0.31, 0.535, 0.6};       /* q-axis open circuit time constant */
48 PetscScalar       M[3];                               /* M = 2*H/w_s */
49 PetscScalar       D[3];                               /* D = 0.1*M */
50 
51 PetscScalar       TM[3]; /* Mechanical Torque */
52 /* Exciter system constants */
53 const PetscScalar KA[3] = {20.0, 20.0, 20.0};    /* Voltage regulartor gain constant */
54 const PetscScalar TA[3] = {0.2, 0.2, 0.2};       /* Voltage regulator time constant */
55 const PetscScalar KE[3] = {1.0, 1.0, 1.0};       /* Exciter gain constant */
56 const PetscScalar TE[3] = {0.314, 0.314, 0.314}; /* Exciter time constant */
57 const PetscScalar KF[3] = {0.063, 0.063, 0.063}; /* Feedback stabilizer gain constant */
58 const PetscScalar TF[3] = {0.35, 0.35, 0.35};    /* Feedback stabilizer time constant */
59 const PetscScalar k1[3] = {0.0039, 0.0039, 0.0039};
60 const PetscScalar k2[3] = {1.555, 1.555, 1.555}; /* k1 and k2 for calculating the saturation function SE = k1*exp(k2*Efd) */
61 
62 PetscScalar       Vref[3];
63 /* Load constants
64   We use a composite load model that describes the load and reactive powers at each time instant as follows
65   P(t) = \sum\limits_{i=0}^ld_nsegsp \ld_alphap_i*P_D0(\frac{V_m(t)}{V_m0})^\ld_betap_i
66   Q(t) = \sum\limits_{i=0}^ld_nsegsq \ld_alphaq_i*Q_D0(\frac{V_m(t)}{V_m0})^\ld_betaq_i
67   where
68     ld_nsegsp,ld_nsegsq - Number of individual load models for real and reactive power loads
69     ld_alphap,ld_alphap - Percentage contribution (weights) or loads
70     P_D0                - Real power load
71     Q_D0                - Reactive power load
72     V_m(t)              - Voltage magnitude at time t
73     V_m0                - Voltage magnitude at t = 0
74     ld_betap, ld_betaq  - exponents describing the load model for real and reactive part
75 
76     Note: All loads have the same characteristic currently.
77 */
78 const PetscScalar PD0[3]       = {1.25, 0.9, 1.0};
79 const PetscScalar QD0[3]       = {0.5, 0.3, 0.35};
80 const PetscInt    ld_nsegsp[3] = {3, 3, 3};
81 const PetscScalar ld_alphap[3] = {1.0, 0.0, 0.0};
82 const PetscScalar ld_betap[3]  = {2.0, 1.0, 0.0};
83 const PetscInt    ld_nsegsq[3] = {3, 3, 3};
84 const PetscScalar ld_alphaq[3] = {1.0, 0.0, 0.0};
85 const PetscScalar ld_betaq[3]  = {2.0, 1.0, 0.0};
86 
87 typedef struct {
88   DM          dmgen, dmnet;        /* DMs to manage generator and network subsystem */
89   DM          dmpgrid;             /* Composite DM to manage the entire power grid */
90   Mat         Ybus;                /* Network admittance matrix */
91   Vec         V0;                  /* Initial voltage vector (Power flow solution) */
92   PetscReal   tfaulton, tfaultoff; /* Fault on and off times */
93   PetscInt    faultbus;            /* Fault bus */
94   PetscScalar Rfault;
95   PetscReal   t0, tmax;
96   PetscInt    neqs_gen, neqs_net, neqs_pgrid;
97   Mat         Sol; /* Matrix to save solution at each time step */
98   PetscInt    stepnum;
99   PetscBool   alg_flg;
100   PetscReal   t;
101   IS          is_diff;        /* indices for differential equations */
102   IS          is_alg;         /* indices for algebraic equations */
103   PetscReal   freq_u, freq_l; /* upper and lower frequency limit */
104   PetscInt    pow;            /* power coefficient used in the cost function */
105   PetscBool   jacp_flg;
106   Mat         J, Jacp;
107   Mat         DRDU, DRDP;
108 } Userctx;
109 
110 /* Converts from machine frame (dq) to network (phase a real,imag) reference frame */
111 PetscErrorCode dq2ri(PetscScalar Fd, PetscScalar Fq, PetscScalar delta, PetscScalar *Fr, PetscScalar *Fi) {
112   PetscFunctionBegin;
113   *Fr = Fd * PetscSinScalar(delta) + Fq * PetscCosScalar(delta);
114   *Fi = -Fd * PetscCosScalar(delta) + Fq * PetscSinScalar(delta);
115   PetscFunctionReturn(0);
116 }
117 
118 /* Converts from network frame ([phase a real,imag) to machine (dq) reference frame */
119 PetscErrorCode ri2dq(PetscScalar Fr, PetscScalar Fi, PetscScalar delta, PetscScalar *Fd, PetscScalar *Fq) {
120   PetscFunctionBegin;
121   *Fd = Fr * PetscSinScalar(delta) - Fi * PetscCosScalar(delta);
122   *Fq = Fr * PetscCosScalar(delta) + Fi * PetscSinScalar(delta);
123   PetscFunctionReturn(0);
124 }
125 
126 /* Saves the solution at each time to a matrix */
127 PetscErrorCode SaveSolution(TS ts) {
128   Userctx           *user;
129   Vec                X;
130   PetscScalar       *mat;
131   const PetscScalar *x;
132   PetscInt           idx;
133   PetscReal          t;
134 
135   PetscFunctionBegin;
136   PetscCall(TSGetApplicationContext(ts, &user));
137   PetscCall(TSGetTime(ts, &t));
138   PetscCall(TSGetSolution(ts, &X));
139   idx = user->stepnum * (user->neqs_pgrid + 1);
140   PetscCall(MatDenseGetArray(user->Sol, &mat));
141   PetscCall(VecGetArrayRead(X, &x));
142   mat[idx] = t;
143   PetscCall(PetscArraycpy(mat + idx + 1, x, user->neqs_pgrid));
144   PetscCall(MatDenseRestoreArray(user->Sol, &mat));
145   PetscCall(VecRestoreArrayRead(X, &x));
146   user->stepnum++;
147   PetscFunctionReturn(0);
148 }
149 
150 PetscErrorCode SetInitialGuess(Vec X, Userctx *user) {
151   Vec          Xgen, Xnet;
152   PetscScalar *xgen, *xnet;
153   PetscInt     i, idx = 0;
154   PetscScalar  Vr, Vi, IGr, IGi, Vm, Vm2;
155   PetscScalar  Eqp, Edp, delta;
156   PetscScalar  Efd, RF, VR; /* Exciter variables */
157   PetscScalar  Id, Iq;      /* Generator dq axis currents */
158   PetscScalar  theta, Vd, Vq, SE;
159 
160   PetscFunctionBegin;
161   M[0] = 2 * H[0] / w_s;
162   M[1] = 2 * H[1] / w_s;
163   M[2] = 2 * H[2] / w_s;
164   D[0] = 0.1 * M[0];
165   D[1] = 0.1 * M[1];
166   D[2] = 0.1 * M[2];
167 
168   PetscCall(DMCompositeGetLocalVectors(user->dmpgrid, &Xgen, &Xnet));
169 
170   /* Network subsystem initialization */
171   PetscCall(VecCopy(user->V0, Xnet));
172 
173   /* Generator subsystem initialization */
174   PetscCall(VecGetArray(Xgen, &xgen));
175   PetscCall(VecGetArray(Xnet, &xnet));
176 
177   for (i = 0; i < ngen; i++) {
178     Vr  = xnet[2 * gbus[i]];     /* Real part of generator terminal voltage */
179     Vi  = xnet[2 * gbus[i] + 1]; /* Imaginary part of the generator terminal voltage */
180     Vm  = PetscSqrtScalar(Vr * Vr + Vi * Vi);
181     Vm2 = Vm * Vm;
182     IGr = (Vr * PG[i] + Vi * QG[i]) / Vm2;
183     IGi = (Vi * PG[i] - Vr * QG[i]) / Vm2;
184 
185     delta = PetscAtan2Real(Vi + Xq[i] * IGr, Vr - Xq[i] * IGi); /* Machine angle */
186 
187     theta = PETSC_PI / 2.0 - delta;
188 
189     Id = IGr * PetscCosScalar(theta) - IGi * PetscSinScalar(theta); /* d-axis stator current */
190     Iq = IGr * PetscSinScalar(theta) + IGi * PetscCosScalar(theta); /* q-axis stator current */
191 
192     Vd = Vr * PetscCosScalar(theta) - Vi * PetscSinScalar(theta);
193     Vq = Vr * PetscSinScalar(theta) + Vi * PetscCosScalar(theta);
194 
195     Edp = Vd + Rs[i] * Id - Xqp[i] * Iq; /* d-axis transient EMF */
196     Eqp = Vq + Rs[i] * Iq + Xdp[i] * Id; /* q-axis transient EMF */
197 
198     TM[i] = PG[i];
199 
200     /* The generator variables are ordered as [Eqp,Edp,delta,w,Id,Iq] */
201     xgen[idx]     = Eqp;
202     xgen[idx + 1] = Edp;
203     xgen[idx + 2] = delta;
204     xgen[idx + 3] = w_s;
205 
206     idx = idx + 4;
207 
208     xgen[idx]     = Id;
209     xgen[idx + 1] = Iq;
210 
211     idx = idx + 2;
212 
213     /* Exciter */
214     Efd = Eqp + (Xd[i] - Xdp[i]) * Id;
215     SE  = k1[i] * PetscExpScalar(k2[i] * Efd);
216     VR  = KE[i] * Efd + SE;
217     RF  = KF[i] * Efd / TF[i];
218 
219     xgen[idx]     = Efd;
220     xgen[idx + 1] = RF;
221     xgen[idx + 2] = VR;
222 
223     Vref[i] = Vm + (VR / KA[i]);
224 
225     idx = idx + 3;
226   }
227 
228   PetscCall(VecRestoreArray(Xgen, &xgen));
229   PetscCall(VecRestoreArray(Xnet, &xnet));
230 
231   /* PetscCall(VecView(Xgen,0)); */
232   PetscCall(DMCompositeGather(user->dmpgrid, INSERT_VALUES, X, Xgen, Xnet));
233   PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid, &Xgen, &Xnet));
234   PetscFunctionReturn(0);
235 }
236 
237 PetscErrorCode InitialGuess(Vec X, Userctx *user, const PetscScalar PGv[]) {
238   Vec          Xgen, Xnet;
239   PetscScalar *xgen, *xnet;
240   PetscInt     i, idx = 0;
241   PetscScalar  Vr, Vi, IGr, IGi, Vm, Vm2;
242   PetscScalar  Eqp, Edp, delta;
243   PetscScalar  Efd, RF, VR; /* Exciter variables */
244   PetscScalar  Id, Iq;      /* Generator dq axis currents */
245   PetscScalar  theta, Vd, Vq, SE;
246 
247   PetscFunctionBegin;
248   M[0] = 2 * H[0] / w_s;
249   M[1] = 2 * H[1] / w_s;
250   M[2] = 2 * H[2] / w_s;
251   D[0] = 0.1 * M[0];
252   D[1] = 0.1 * M[1];
253   D[2] = 0.1 * M[2];
254 
255   PetscCall(DMCompositeGetLocalVectors(user->dmpgrid, &Xgen, &Xnet));
256 
257   /* Network subsystem initialization */
258   PetscCall(VecCopy(user->V0, Xnet));
259 
260   /* Generator subsystem initialization */
261   PetscCall(VecGetArray(Xgen, &xgen));
262   PetscCall(VecGetArray(Xnet, &xnet));
263 
264   for (i = 0; i < ngen; i++) {
265     Vr  = xnet[2 * gbus[i]];     /* Real part of generator terminal voltage */
266     Vi  = xnet[2 * gbus[i] + 1]; /* Imaginary part of the generator terminal voltage */
267     Vm  = PetscSqrtScalar(Vr * Vr + Vi * Vi);
268     Vm2 = Vm * Vm;
269     IGr = (Vr * PGv[i] + Vi * QG[i]) / Vm2;
270     IGi = (Vi * PGv[i] - Vr * QG[i]) / Vm2;
271 
272     delta = PetscAtan2Real(Vi + Xq[i] * IGr, Vr - Xq[i] * IGi); /* Machine angle */
273 
274     theta = PETSC_PI / 2.0 - delta;
275 
276     Id = IGr * PetscCosScalar(theta) - IGi * PetscSinScalar(theta); /* d-axis stator current */
277     Iq = IGr * PetscSinScalar(theta) + IGi * PetscCosScalar(theta); /* q-axis stator current */
278 
279     Vd = Vr * PetscCosScalar(theta) - Vi * PetscSinScalar(theta);
280     Vq = Vr * PetscSinScalar(theta) + Vi * PetscCosScalar(theta);
281 
282     Edp = Vd + Rs[i] * Id - Xqp[i] * Iq; /* d-axis transient EMF */
283     Eqp = Vq + Rs[i] * Iq + Xdp[i] * Id; /* q-axis transient EMF */
284 
285     /* The generator variables are ordered as [Eqp,Edp,delta,w,Id,Iq] */
286     xgen[idx]     = Eqp;
287     xgen[idx + 1] = Edp;
288     xgen[idx + 2] = delta;
289     xgen[idx + 3] = w_s;
290 
291     idx = idx + 4;
292 
293     xgen[idx]     = Id;
294     xgen[idx + 1] = Iq;
295 
296     idx = idx + 2;
297 
298     /* Exciter */
299     Efd = Eqp + (Xd[i] - Xdp[i]) * Id;
300     SE  = k1[i] * PetscExpScalar(k2[i] * Efd);
301     VR  = KE[i] * Efd + SE;
302     RF  = KF[i] * Efd / TF[i];
303 
304     xgen[idx]     = Efd;
305     xgen[idx + 1] = RF;
306     xgen[idx + 2] = VR;
307 
308     idx = idx + 3;
309   }
310 
311   PetscCall(VecRestoreArray(Xgen, &xgen));
312   PetscCall(VecRestoreArray(Xnet, &xnet));
313 
314   /* PetscCall(VecView(Xgen,0)); */
315   PetscCall(DMCompositeGather(user->dmpgrid, INSERT_VALUES, X, Xgen, Xnet));
316   PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid, &Xgen, &Xnet));
317   PetscFunctionReturn(0);
318 }
319 
320 PetscErrorCode DICDPFiniteDifference(Vec X, Vec *DICDP, Userctx *user) {
321   Vec         Y;
322   PetscScalar PGv[3], eps;
323   PetscInt    i, j;
324 
325   PetscFunctionBegin;
326   eps = 1.e-7;
327   PetscCall(VecDuplicate(X, &Y));
328 
329   for (i = 0; i < ngen; i++) {
330     for (j = 0; j < 3; j++) PGv[j] = PG[j];
331     PGv[i] = PG[i] + eps;
332     PetscCall(InitialGuess(Y, user, PGv));
333     PetscCall(InitialGuess(X, user, PG));
334 
335     PetscCall(VecAXPY(Y, -1.0, X));
336     PetscCall(VecScale(Y, 1. / eps));
337     PetscCall(VecCopy(Y, DICDP[i]));
338   }
339   PetscCall(VecDestroy(&Y));
340   PetscFunctionReturn(0);
341 }
342 
343 /* Computes F = [-f(x,y);g(x,y)] */
344 PetscErrorCode ResidualFunction(SNES snes, Vec X, Vec F, Userctx *user) {
345   Vec          Xgen, Xnet, Fgen, Fnet;
346   PetscScalar *xgen, *xnet, *fgen, *fnet;
347   PetscInt     i, idx = 0;
348   PetscScalar  Vr, Vi, Vm, Vm2;
349   PetscScalar  Eqp, Edp, delta, w; /* Generator variables */
350   PetscScalar  Efd, RF, VR;        /* Exciter variables */
351   PetscScalar  Id, Iq;             /* Generator dq axis currents */
352   PetscScalar  Vd, Vq, SE;
353   PetscScalar  IGr, IGi, IDr, IDi;
354   PetscScalar  Zdq_inv[4], det;
355   PetscScalar  PD, QD, Vm0, *v0;
356   PetscInt     k;
357 
358   PetscFunctionBegin;
359   PetscCall(VecZeroEntries(F));
360   PetscCall(DMCompositeGetLocalVectors(user->dmpgrid, &Xgen, &Xnet));
361   PetscCall(DMCompositeGetLocalVectors(user->dmpgrid, &Fgen, &Fnet));
362   PetscCall(DMCompositeScatter(user->dmpgrid, X, Xgen, Xnet));
363   PetscCall(DMCompositeScatter(user->dmpgrid, F, Fgen, Fnet));
364 
365   /* Network current balance residual IG + Y*V + IL = 0. Only YV is added here.
366      The generator current injection, IG, and load current injection, ID are added later
367   */
368   /* Note that the values in Ybus are stored assuming the imaginary current balance
369      equation is ordered first followed by real current balance equation for each bus.
370      Thus imaginary current contribution goes in location 2*i, and
371      real current contribution in 2*i+1
372   */
373   PetscCall(MatMult(user->Ybus, Xnet, Fnet));
374 
375   PetscCall(VecGetArray(Xgen, &xgen));
376   PetscCall(VecGetArray(Xnet, &xnet));
377   PetscCall(VecGetArray(Fgen, &fgen));
378   PetscCall(VecGetArray(Fnet, &fnet));
379 
380   /* Generator subsystem */
381   for (i = 0; i < ngen; i++) {
382     Eqp   = xgen[idx];
383     Edp   = xgen[idx + 1];
384     delta = xgen[idx + 2];
385     w     = xgen[idx + 3];
386     Id    = xgen[idx + 4];
387     Iq    = xgen[idx + 5];
388     Efd   = xgen[idx + 6];
389     RF    = xgen[idx + 7];
390     VR    = xgen[idx + 8];
391 
392     /* Generator differential equations */
393     fgen[idx]     = (Eqp + (Xd[i] - Xdp[i]) * Id - Efd) / Td0p[i];
394     fgen[idx + 1] = (Edp - (Xq[i] - Xqp[i]) * Iq) / Tq0p[i];
395     fgen[idx + 2] = -w + w_s;
396     fgen[idx + 3] = (-TM[i] + Edp * Id + Eqp * Iq + (Xqp[i] - Xdp[i]) * Id * Iq + D[i] * (w - w_s)) / M[i];
397 
398     Vr = xnet[2 * gbus[i]];     /* Real part of generator terminal voltage */
399     Vi = xnet[2 * gbus[i] + 1]; /* Imaginary part of the generator terminal voltage */
400 
401     PetscCall(ri2dq(Vr, Vi, delta, &Vd, &Vq));
402     /* Algebraic equations for stator currents */
403     det = Rs[i] * Rs[i] + Xdp[i] * Xqp[i];
404 
405     Zdq_inv[0] = Rs[i] / det;
406     Zdq_inv[1] = Xqp[i] / det;
407     Zdq_inv[2] = -Xdp[i] / det;
408     Zdq_inv[3] = Rs[i] / det;
409 
410     fgen[idx + 4] = Zdq_inv[0] * (-Edp + Vd) + Zdq_inv[1] * (-Eqp + Vq) + Id;
411     fgen[idx + 5] = Zdq_inv[2] * (-Edp + Vd) + Zdq_inv[3] * (-Eqp + Vq) + Iq;
412 
413     /* Add generator current injection to network */
414     PetscCall(dq2ri(Id, Iq, delta, &IGr, &IGi));
415 
416     fnet[2 * gbus[i]] -= IGi;
417     fnet[2 * gbus[i] + 1] -= IGr;
418 
419     Vm = PetscSqrtScalar(Vd * Vd + Vq * Vq);
420 
421     SE = k1[i] * PetscExpScalar(k2[i] * Efd);
422 
423     /* Exciter differential equations */
424     fgen[idx + 6] = (KE[i] * Efd + SE - VR) / TE[i];
425     fgen[idx + 7] = (RF - KF[i] * Efd / TF[i]) / TF[i];
426     fgen[idx + 8] = (VR - KA[i] * RF + KA[i] * KF[i] * Efd / TF[i] - KA[i] * (Vref[i] - Vm)) / TA[i];
427 
428     idx = idx + 9;
429   }
430 
431   PetscCall(VecGetArray(user->V0, &v0));
432   for (i = 0; i < nload; i++) {
433     Vr  = xnet[2 * lbus[i]];     /* Real part of load bus voltage */
434     Vi  = xnet[2 * lbus[i] + 1]; /* Imaginary part of the load bus voltage */
435     Vm  = PetscSqrtScalar(Vr * Vr + Vi * Vi);
436     Vm2 = Vm * Vm;
437     Vm0 = PetscSqrtScalar(v0[2 * lbus[i]] * v0[2 * lbus[i]] + v0[2 * lbus[i] + 1] * v0[2 * lbus[i] + 1]);
438     PD = QD = 0.0;
439     for (k = 0; k < ld_nsegsp[i]; k++) PD += ld_alphap[k] * PD0[i] * PetscPowScalar((Vm / Vm0), ld_betap[k]);
440     for (k = 0; k < ld_nsegsq[i]; k++) QD += ld_alphaq[k] * QD0[i] * PetscPowScalar((Vm / Vm0), ld_betaq[k]);
441 
442     /* Load currents */
443     IDr = (PD * Vr + QD * Vi) / Vm2;
444     IDi = (-QD * Vr + PD * Vi) / Vm2;
445 
446     fnet[2 * lbus[i]] += IDi;
447     fnet[2 * lbus[i] + 1] += IDr;
448   }
449   PetscCall(VecRestoreArray(user->V0, &v0));
450 
451   PetscCall(VecRestoreArray(Xgen, &xgen));
452   PetscCall(VecRestoreArray(Xnet, &xnet));
453   PetscCall(VecRestoreArray(Fgen, &fgen));
454   PetscCall(VecRestoreArray(Fnet, &fnet));
455 
456   PetscCall(DMCompositeGather(user->dmpgrid, INSERT_VALUES, F, Fgen, Fnet));
457   PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid, &Xgen, &Xnet));
458   PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid, &Fgen, &Fnet));
459   PetscFunctionReturn(0);
460 }
461 
462 /* \dot{x} - f(x,y)
463      g(x,y) = 0
464  */
465 PetscErrorCode IFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, Userctx *user) {
466   SNES               snes;
467   PetscScalar       *f;
468   const PetscScalar *xdot;
469   PetscInt           i;
470 
471   PetscFunctionBegin;
472   user->t = t;
473 
474   PetscCall(TSGetSNES(ts, &snes));
475   PetscCall(ResidualFunction(snes, X, F, user));
476   PetscCall(VecGetArray(F, &f));
477   PetscCall(VecGetArrayRead(Xdot, &xdot));
478   for (i = 0; i < ngen; i++) {
479     f[9 * i] += xdot[9 * i];
480     f[9 * i + 1] += xdot[9 * i + 1];
481     f[9 * i + 2] += xdot[9 * i + 2];
482     f[9 * i + 3] += xdot[9 * i + 3];
483     f[9 * i + 6] += xdot[9 * i + 6];
484     f[9 * i + 7] += xdot[9 * i + 7];
485     f[9 * i + 8] += xdot[9 * i + 8];
486   }
487   PetscCall(VecRestoreArray(F, &f));
488   PetscCall(VecRestoreArrayRead(Xdot, &xdot));
489   PetscFunctionReturn(0);
490 }
491 
492 /* This function is used for solving the algebraic system only during fault on and
493    off times. It computes the entire F and then zeros out the part corresponding to
494    differential equations
495  F = [0;g(y)];
496 */
497 PetscErrorCode AlgFunction(SNES snes, Vec X, Vec F, void *ctx) {
498   Userctx     *user = (Userctx *)ctx;
499   PetscScalar *f;
500   PetscInt     i;
501 
502   PetscFunctionBegin;
503   PetscCall(ResidualFunction(snes, X, F, user));
504   PetscCall(VecGetArray(F, &f));
505   for (i = 0; i < ngen; i++) {
506     f[9 * i]     = 0;
507     f[9 * i + 1] = 0;
508     f[9 * i + 2] = 0;
509     f[9 * i + 3] = 0;
510     f[9 * i + 6] = 0;
511     f[9 * i + 7] = 0;
512     f[9 * i + 8] = 0;
513   }
514   PetscCall(VecRestoreArray(F, &f));
515   PetscFunctionReturn(0);
516 }
517 
518 PetscErrorCode PreallocateJacobian(Mat J, Userctx *user) {
519   PetscInt *d_nnz;
520   PetscInt  i, idx = 0, start = 0;
521   PetscInt  ncols;
522 
523   PetscFunctionBegin;
524   PetscCall(PetscMalloc1(user->neqs_pgrid, &d_nnz));
525   for (i = 0; i < user->neqs_pgrid; i++) d_nnz[i] = 0;
526   /* Generator subsystem */
527   for (i = 0; i < ngen; i++) {
528     d_nnz[idx] += 3;
529     d_nnz[idx + 1] += 2;
530     d_nnz[idx + 2] += 2;
531     d_nnz[idx + 3] += 5;
532     d_nnz[idx + 4] += 6;
533     d_nnz[idx + 5] += 6;
534 
535     d_nnz[user->neqs_gen + 2 * gbus[i]] += 3;
536     d_nnz[user->neqs_gen + 2 * gbus[i] + 1] += 3;
537 
538     d_nnz[idx + 6] += 2;
539     d_nnz[idx + 7] += 2;
540     d_nnz[idx + 8] += 5;
541 
542     idx = idx + 9;
543   }
544 
545   start = user->neqs_gen;
546   for (i = 0; i < nbus; i++) {
547     PetscCall(MatGetRow(user->Ybus, 2 * i, &ncols, NULL, NULL));
548     d_nnz[start + 2 * i] += ncols;
549     d_nnz[start + 2 * i + 1] += ncols;
550     PetscCall(MatRestoreRow(user->Ybus, 2 * i, &ncols, NULL, NULL));
551   }
552 
553   PetscCall(MatSeqAIJSetPreallocation(J, 0, d_nnz));
554   PetscCall(PetscFree(d_nnz));
555   PetscFunctionReturn(0);
556 }
557 
558 /*
559    J = [-df_dx, -df_dy
560         dg_dx, dg_dy]
561 */
562 PetscErrorCode ResidualJacobian(SNES snes, Vec X, Mat J, Mat B, void *ctx) {
563   Userctx           *user = (Userctx *)ctx;
564   Vec                Xgen, Xnet;
565   PetscScalar       *xgen, *xnet;
566   PetscInt           i, idx = 0;
567   PetscScalar        Vr, Vi, Vm, Vm2;
568   PetscScalar        Eqp, Edp, delta; /* Generator variables */
569   PetscScalar        Efd;             /* Exciter variables */
570   PetscScalar        Id, Iq;          /* Generator dq axis currents */
571   PetscScalar        Vd, Vq;
572   PetscScalar        val[10];
573   PetscInt           row[2], col[10];
574   PetscInt           net_start = user->neqs_gen;
575   PetscInt           ncols;
576   const PetscInt    *cols;
577   const PetscScalar *yvals;
578   PetscInt           k;
579   PetscScalar        Zdq_inv[4], det;
580   PetscScalar        dVd_dVr, dVd_dVi, dVq_dVr, dVq_dVi, dVd_ddelta, dVq_ddelta;
581   PetscScalar        dIGr_ddelta, dIGi_ddelta, dIGr_dId, dIGr_dIq, dIGi_dId, dIGi_dIq;
582   PetscScalar        dSE_dEfd;
583   PetscScalar        dVm_dVd, dVm_dVq, dVm_dVr, dVm_dVi;
584   PetscScalar        PD, QD, Vm0, *v0, Vm4;
585   PetscScalar        dPD_dVr, dPD_dVi, dQD_dVr, dQD_dVi;
586   PetscScalar        dIDr_dVr, dIDr_dVi, dIDi_dVr, dIDi_dVi;
587 
588   PetscFunctionBegin;
589   PetscCall(MatZeroEntries(B));
590   PetscCall(DMCompositeGetLocalVectors(user->dmpgrid, &Xgen, &Xnet));
591   PetscCall(DMCompositeScatter(user->dmpgrid, X, Xgen, Xnet));
592 
593   PetscCall(VecGetArray(Xgen, &xgen));
594   PetscCall(VecGetArray(Xnet, &xnet));
595 
596   /* Generator subsystem */
597   for (i = 0; i < ngen; i++) {
598     Eqp   = xgen[idx];
599     Edp   = xgen[idx + 1];
600     delta = xgen[idx + 2];
601     Id    = xgen[idx + 4];
602     Iq    = xgen[idx + 5];
603     Efd   = xgen[idx + 6];
604 
605     /*    fgen[idx]   = (Eqp + (Xd[i] - Xdp[i])*Id - Efd)/Td0p[i]; */
606     row[0] = idx;
607     col[0] = idx;
608     col[1] = idx + 4;
609     col[2] = idx + 6;
610     val[0] = 1 / Td0p[i];
611     val[1] = (Xd[i] - Xdp[i]) / Td0p[i];
612     val[2] = -1 / Td0p[i];
613 
614     PetscCall(MatSetValues(J, 1, row, 3, col, val, INSERT_VALUES));
615 
616     /*    fgen[idx+1] = (Edp - (Xq[i] - Xqp[i])*Iq)/Tq0p[i]; */
617     row[0] = idx + 1;
618     col[0] = idx + 1;
619     col[1] = idx + 5;
620     val[0] = 1 / Tq0p[i];
621     val[1] = -(Xq[i] - Xqp[i]) / Tq0p[i];
622     PetscCall(MatSetValues(J, 1, row, 2, col, val, INSERT_VALUES));
623 
624     /*    fgen[idx+2] = - w + w_s; */
625     row[0] = idx + 2;
626     col[0] = idx + 2;
627     col[1] = idx + 3;
628     val[0] = 0;
629     val[1] = -1;
630     PetscCall(MatSetValues(J, 1, row, 2, col, val, INSERT_VALUES));
631 
632     /*    fgen[idx+3] = (-TM[i] + Edp*Id + Eqp*Iq + (Xqp[i] - Xdp[i])*Id*Iq + D[i]*(w - w_s))/M[i]; */
633     row[0] = idx + 3;
634     col[0] = idx;
635     col[1] = idx + 1;
636     col[2] = idx + 3;
637     col[3] = idx + 4;
638     col[4] = idx + 5;
639     val[0] = Iq / M[i];
640     val[1] = Id / M[i];
641     val[2] = D[i] / M[i];
642     val[3] = (Edp + (Xqp[i] - Xdp[i]) * Iq) / M[i];
643     val[4] = (Eqp + (Xqp[i] - Xdp[i]) * Id) / M[i];
644     PetscCall(MatSetValues(J, 1, row, 5, col, val, INSERT_VALUES));
645 
646     Vr = xnet[2 * gbus[i]];     /* Real part of generator terminal voltage */
647     Vi = xnet[2 * gbus[i] + 1]; /* Imaginary part of the generator terminal voltage */
648     PetscCall(ri2dq(Vr, Vi, delta, &Vd, &Vq));
649 
650     det = Rs[i] * Rs[i] + Xdp[i] * Xqp[i];
651 
652     Zdq_inv[0] = Rs[i] / det;
653     Zdq_inv[1] = Xqp[i] / det;
654     Zdq_inv[2] = -Xdp[i] / det;
655     Zdq_inv[3] = Rs[i] / det;
656 
657     dVd_dVr    = PetscSinScalar(delta);
658     dVd_dVi    = -PetscCosScalar(delta);
659     dVq_dVr    = PetscCosScalar(delta);
660     dVq_dVi    = PetscSinScalar(delta);
661     dVd_ddelta = Vr * PetscCosScalar(delta) + Vi * PetscSinScalar(delta);
662     dVq_ddelta = -Vr * PetscSinScalar(delta) + Vi * PetscCosScalar(delta);
663 
664     /*    fgen[idx+4] = Zdq_inv[0]*(-Edp + Vd) + Zdq_inv[1]*(-Eqp + Vq) + Id; */
665     row[0] = idx + 4;
666     col[0] = idx;
667     col[1] = idx + 1;
668     col[2] = idx + 2;
669     val[0] = -Zdq_inv[1];
670     val[1] = -Zdq_inv[0];
671     val[2] = Zdq_inv[0] * dVd_ddelta + Zdq_inv[1] * dVq_ddelta;
672     col[3] = idx + 4;
673     col[4] = net_start + 2 * gbus[i];
674     col[5] = net_start + 2 * gbus[i] + 1;
675     val[3] = 1;
676     val[4] = Zdq_inv[0] * dVd_dVr + Zdq_inv[1] * dVq_dVr;
677     val[5] = Zdq_inv[0] * dVd_dVi + Zdq_inv[1] * dVq_dVi;
678     PetscCall(MatSetValues(J, 1, row, 6, col, val, INSERT_VALUES));
679 
680     /*  fgen[idx+5] = Zdq_inv[2]*(-Edp + Vd) + Zdq_inv[3]*(-Eqp + Vq) + Iq; */
681     row[0] = idx + 5;
682     col[0] = idx;
683     col[1] = idx + 1;
684     col[2] = idx + 2;
685     val[0] = -Zdq_inv[3];
686     val[1] = -Zdq_inv[2];
687     val[2] = Zdq_inv[2] * dVd_ddelta + Zdq_inv[3] * dVq_ddelta;
688     col[3] = idx + 5;
689     col[4] = net_start + 2 * gbus[i];
690     col[5] = net_start + 2 * gbus[i] + 1;
691     val[3] = 1;
692     val[4] = Zdq_inv[2] * dVd_dVr + Zdq_inv[3] * dVq_dVr;
693     val[5] = Zdq_inv[2] * dVd_dVi + Zdq_inv[3] * dVq_dVi;
694     PetscCall(MatSetValues(J, 1, row, 6, col, val, INSERT_VALUES));
695 
696     dIGr_ddelta = Id * PetscCosScalar(delta) - Iq * PetscSinScalar(delta);
697     dIGi_ddelta = Id * PetscSinScalar(delta) + Iq * PetscCosScalar(delta);
698     dIGr_dId    = PetscSinScalar(delta);
699     dIGr_dIq    = PetscCosScalar(delta);
700     dIGi_dId    = -PetscCosScalar(delta);
701     dIGi_dIq    = PetscSinScalar(delta);
702 
703     /* fnet[2*gbus[i]]   -= IGi; */
704     row[0] = net_start + 2 * gbus[i];
705     col[0] = idx + 2;
706     col[1] = idx + 4;
707     col[2] = idx + 5;
708     val[0] = -dIGi_ddelta;
709     val[1] = -dIGi_dId;
710     val[2] = -dIGi_dIq;
711     PetscCall(MatSetValues(J, 1, row, 3, col, val, INSERT_VALUES));
712 
713     /* fnet[2*gbus[i]+1]   -= IGr; */
714     row[0] = net_start + 2 * gbus[i] + 1;
715     col[0] = idx + 2;
716     col[1] = idx + 4;
717     col[2] = idx + 5;
718     val[0] = -dIGr_ddelta;
719     val[1] = -dIGr_dId;
720     val[2] = -dIGr_dIq;
721     PetscCall(MatSetValues(J, 1, row, 3, col, val, INSERT_VALUES));
722 
723     Vm = PetscSqrtScalar(Vd * Vd + Vq * Vq);
724 
725     /*    fgen[idx+6] = (KE[i]*Efd + SE - VR)/TE[i]; */
726     /*    SE  = k1[i]*PetscExpScalar(k2[i]*Efd); */
727     dSE_dEfd = k1[i] * k2[i] * PetscExpScalar(k2[i] * Efd);
728 
729     row[0] = idx + 6;
730     col[0] = idx + 6;
731     col[1] = idx + 8;
732     val[0] = (KE[i] + dSE_dEfd) / TE[i];
733     val[1] = -1 / TE[i];
734     PetscCall(MatSetValues(J, 1, row, 2, col, val, INSERT_VALUES));
735 
736     /* Exciter differential equations */
737 
738     /*    fgen[idx+7] = (RF - KF[i]*Efd/TF[i])/TF[i]; */
739     row[0] = idx + 7;
740     col[0] = idx + 6;
741     col[1] = idx + 7;
742     val[0] = (-KF[i] / TF[i]) / TF[i];
743     val[1] = 1 / TF[i];
744     PetscCall(MatSetValues(J, 1, row, 2, col, val, INSERT_VALUES));
745 
746     /*    fgen[idx+8] = (VR - KA[i]*RF + KA[i]*KF[i]*Efd/TF[i] - KA[i]*(Vref[i] - Vm))/TA[i]; */
747     /* Vm = (Vd^2 + Vq^2)^0.5; */
748     dVm_dVd = Vd / Vm;
749     dVm_dVq = Vq / Vm;
750     dVm_dVr = dVm_dVd * dVd_dVr + dVm_dVq * dVq_dVr;
751     dVm_dVi = dVm_dVd * dVd_dVi + dVm_dVq * dVq_dVi;
752     row[0]  = idx + 8;
753     col[0]  = idx + 6;
754     col[1]  = idx + 7;
755     col[2]  = idx + 8;
756     val[0]  = (KA[i] * KF[i] / TF[i]) / TA[i];
757     val[1]  = -KA[i] / TA[i];
758     val[2]  = 1 / TA[i];
759     col[3]  = net_start + 2 * gbus[i];
760     col[4]  = net_start + 2 * gbus[i] + 1;
761     val[3]  = KA[i] * dVm_dVr / TA[i];
762     val[4]  = KA[i] * dVm_dVi / TA[i];
763     PetscCall(MatSetValues(J, 1, row, 5, col, val, INSERT_VALUES));
764     idx = idx + 9;
765   }
766 
767   for (i = 0; i < nbus; i++) {
768     PetscCall(MatGetRow(user->Ybus, 2 * i, &ncols, &cols, &yvals));
769     row[0] = net_start + 2 * i;
770     for (k = 0; k < ncols; k++) {
771       col[k] = net_start + cols[k];
772       val[k] = yvals[k];
773     }
774     PetscCall(MatSetValues(J, 1, row, ncols, col, val, INSERT_VALUES));
775     PetscCall(MatRestoreRow(user->Ybus, 2 * i, &ncols, &cols, &yvals));
776 
777     PetscCall(MatGetRow(user->Ybus, 2 * i + 1, &ncols, &cols, &yvals));
778     row[0] = net_start + 2 * i + 1;
779     for (k = 0; k < ncols; k++) {
780       col[k] = net_start + cols[k];
781       val[k] = yvals[k];
782     }
783     PetscCall(MatSetValues(J, 1, row, ncols, col, val, INSERT_VALUES));
784     PetscCall(MatRestoreRow(user->Ybus, 2 * i + 1, &ncols, &cols, &yvals));
785   }
786 
787   PetscCall(MatAssemblyBegin(J, MAT_FLUSH_ASSEMBLY));
788   PetscCall(MatAssemblyEnd(J, MAT_FLUSH_ASSEMBLY));
789 
790   PetscCall(VecGetArray(user->V0, &v0));
791   for (i = 0; i < nload; i++) {
792     Vr  = xnet[2 * lbus[i]];     /* Real part of load bus voltage */
793     Vi  = xnet[2 * lbus[i] + 1]; /* Imaginary part of the load bus voltage */
794     Vm  = PetscSqrtScalar(Vr * Vr + Vi * Vi);
795     Vm2 = Vm * Vm;
796     Vm4 = Vm2 * Vm2;
797     Vm0 = PetscSqrtScalar(v0[2 * lbus[i]] * v0[2 * lbus[i]] + v0[2 * lbus[i] + 1] * v0[2 * lbus[i] + 1]);
798     PD = QD = 0.0;
799     dPD_dVr = dPD_dVi = dQD_dVr = dQD_dVi = 0.0;
800     for (k = 0; k < ld_nsegsp[i]; k++) {
801       PD += ld_alphap[k] * PD0[i] * PetscPowScalar((Vm / Vm0), ld_betap[k]);
802       dPD_dVr += ld_alphap[k] * ld_betap[k] * PD0[i] * PetscPowScalar((1 / Vm0), ld_betap[k]) * Vr * PetscPowScalar(Vm, (ld_betap[k] - 2));
803       dPD_dVi += ld_alphap[k] * ld_betap[k] * PD0[i] * PetscPowScalar((1 / Vm0), ld_betap[k]) * Vi * PetscPowScalar(Vm, (ld_betap[k] - 2));
804     }
805     for (k = 0; k < ld_nsegsq[i]; k++) {
806       QD += ld_alphaq[k] * QD0[i] * PetscPowScalar((Vm / Vm0), ld_betaq[k]);
807       dQD_dVr += ld_alphaq[k] * ld_betaq[k] * QD0[i] * PetscPowScalar((1 / Vm0), ld_betaq[k]) * Vr * PetscPowScalar(Vm, (ld_betaq[k] - 2));
808       dQD_dVi += ld_alphaq[k] * ld_betaq[k] * QD0[i] * PetscPowScalar((1 / Vm0), ld_betaq[k]) * Vi * PetscPowScalar(Vm, (ld_betaq[k] - 2));
809     }
810 
811     /*    IDr = (PD*Vr + QD*Vi)/Vm2; */
812     /*    IDi = (-QD*Vr + PD*Vi)/Vm2; */
813 
814     dIDr_dVr = (dPD_dVr * Vr + dQD_dVr * Vi + PD) / Vm2 - ((PD * Vr + QD * Vi) * 2 * Vr) / Vm4;
815     dIDr_dVi = (dPD_dVi * Vr + dQD_dVi * Vi + QD) / Vm2 - ((PD * Vr + QD * Vi) * 2 * Vi) / Vm4;
816 
817     dIDi_dVr = (-dQD_dVr * Vr + dPD_dVr * Vi - QD) / Vm2 - ((-QD * Vr + PD * Vi) * 2 * Vr) / Vm4;
818     dIDi_dVi = (-dQD_dVi * Vr + dPD_dVi * Vi + PD) / Vm2 - ((-QD * Vr + PD * Vi) * 2 * Vi) / Vm4;
819 
820     /*    fnet[2*lbus[i]]   += IDi; */
821     row[0] = net_start + 2 * lbus[i];
822     col[0] = net_start + 2 * lbus[i];
823     col[1] = net_start + 2 * lbus[i] + 1;
824     val[0] = dIDi_dVr;
825     val[1] = dIDi_dVi;
826     PetscCall(MatSetValues(J, 1, row, 2, col, val, ADD_VALUES));
827     /*    fnet[2*lbus[i]+1] += IDr; */
828     row[0] = net_start + 2 * lbus[i] + 1;
829     col[0] = net_start + 2 * lbus[i];
830     col[1] = net_start + 2 * lbus[i] + 1;
831     val[0] = dIDr_dVr;
832     val[1] = dIDr_dVi;
833     PetscCall(MatSetValues(J, 1, row, 2, col, val, ADD_VALUES));
834   }
835   PetscCall(VecRestoreArray(user->V0, &v0));
836 
837   PetscCall(VecRestoreArray(Xgen, &xgen));
838   PetscCall(VecRestoreArray(Xnet, &xnet));
839 
840   PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid, &Xgen, &Xnet));
841 
842   PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
843   PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
844   PetscFunctionReturn(0);
845 }
846 
847 /*
848    J = [I, 0
849         dg_dx, dg_dy]
850 */
851 PetscErrorCode AlgJacobian(SNES snes, Vec X, Mat A, Mat B, void *ctx) {
852   Userctx *user = (Userctx *)ctx;
853 
854   PetscFunctionBegin;
855   PetscCall(ResidualJacobian(snes, X, A, B, ctx));
856   PetscCall(MatSetOption(A, MAT_KEEP_NONZERO_PATTERN, PETSC_TRUE));
857   PetscCall(MatZeroRowsIS(A, user->is_diff, 1.0, NULL, NULL));
858   PetscFunctionReturn(0);
859 }
860 
861 /*
862    J = [a*I-df_dx, -df_dy
863         dg_dx, dg_dy]
864 */
865 
866 PetscErrorCode IJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat A, Mat B, Userctx *user) {
867   SNES        snes;
868   PetscScalar atmp = (PetscScalar)a;
869   PetscInt    i, row;
870 
871   PetscFunctionBegin;
872   user->t = t;
873 
874   PetscCall(TSGetSNES(ts, &snes));
875   PetscCall(ResidualJacobian(snes, X, A, B, user));
876   for (i = 0; i < ngen; i++) {
877     row = 9 * i;
878     PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES));
879     row = 9 * i + 1;
880     PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES));
881     row = 9 * i + 2;
882     PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES));
883     row = 9 * i + 3;
884     PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES));
885     row = 9 * i + 6;
886     PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES));
887     row = 9 * i + 7;
888     PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES));
889     row = 9 * i + 8;
890     PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES));
891   }
892   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
893   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
894   PetscFunctionReturn(0);
895 }
896 
897 /* Matrix JacobianP is constant so that it only needs to be evaluated once */
898 static PetscErrorCode RHSJacobianP(TS ts, PetscReal t, Vec X, Mat A, void *ctx0) {
899   PetscScalar a;
900   PetscInt    row, col;
901   Userctx    *ctx = (Userctx *)ctx0;
902 
903   PetscFunctionBeginUser;
904 
905   if (ctx->jacp_flg) {
906     PetscCall(MatZeroEntries(A));
907 
908     for (col = 0; col < 3; col++) {
909       a   = 1.0 / M[col];
910       row = 9 * col + 3;
911       PetscCall(MatSetValues(A, 1, &row, 1, &col, &a, INSERT_VALUES));
912     }
913 
914     ctx->jacp_flg = PETSC_FALSE;
915 
916     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
917     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
918   }
919   PetscFunctionReturn(0);
920 }
921 
922 static PetscErrorCode CostIntegrand(TS ts, PetscReal t, Vec U, Vec R, Userctx *user) {
923   const PetscScalar *u;
924   PetscInt           idx;
925   Vec                Xgen, Xnet;
926   PetscScalar       *r, *xgen;
927   PetscInt           i;
928 
929   PetscFunctionBegin;
930   PetscCall(DMCompositeGetLocalVectors(user->dmpgrid, &Xgen, &Xnet));
931   PetscCall(DMCompositeScatter(user->dmpgrid, U, Xgen, Xnet));
932 
933   PetscCall(VecGetArray(Xgen, &xgen));
934 
935   PetscCall(VecGetArrayRead(U, &u));
936   PetscCall(VecGetArray(R, &r));
937   r[0] = 0.;
938   idx  = 0;
939   for (i = 0; i < ngen; i++) {
940     r[0] += PetscPowScalarInt(PetscMax(0., PetscMax(xgen[idx + 3] / (2. * PETSC_PI) - user->freq_u, user->freq_l - xgen[idx + 3] / (2. * PETSC_PI))), user->pow);
941     idx += 9;
942   }
943   PetscCall(VecRestoreArrayRead(U, &u));
944   PetscCall(VecRestoreArray(R, &r));
945   PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid, &Xgen, &Xnet));
946   PetscFunctionReturn(0);
947 }
948 
949 static PetscErrorCode DRDUJacobianTranspose(TS ts, PetscReal t, Vec U, Mat DRDU, Mat B, Userctx *user) {
950   Vec          Xgen, Xnet, Dgen, Dnet;
951   PetscScalar *xgen, *dgen;
952   PetscInt     i;
953   PetscInt     idx;
954   Vec          drdu_col;
955   PetscScalar *xarr;
956 
957   PetscFunctionBegin;
958   PetscCall(VecDuplicate(U, &drdu_col));
959   PetscCall(MatDenseGetColumn(DRDU, 0, &xarr));
960   PetscCall(VecPlaceArray(drdu_col, xarr));
961   PetscCall(DMCompositeGetLocalVectors(user->dmpgrid, &Xgen, &Xnet));
962   PetscCall(DMCompositeGetLocalVectors(user->dmpgrid, &Dgen, &Dnet));
963   PetscCall(DMCompositeScatter(user->dmpgrid, U, Xgen, Xnet));
964   PetscCall(DMCompositeScatter(user->dmpgrid, drdu_col, Dgen, Dnet));
965 
966   PetscCall(VecGetArray(Xgen, &xgen));
967   PetscCall(VecGetArray(Dgen, &dgen));
968 
969   idx = 0;
970   for (i = 0; i < ngen; i++) {
971     dgen[idx + 3] = 0.;
972     if (xgen[idx + 3] / (2. * PETSC_PI) > user->freq_u) dgen[idx + 3] = user->pow * PetscPowScalarInt(xgen[idx + 3] / (2. * PETSC_PI) - user->freq_u, user->pow - 1) / (2. * PETSC_PI);
973     if (xgen[idx + 3] / (2. * PETSC_PI) < user->freq_l) dgen[idx + 3] = user->pow * PetscPowScalarInt(user->freq_l - xgen[idx + 3] / (2. * PETSC_PI), user->pow - 1) / (-2. * PETSC_PI);
974     idx += 9;
975   }
976 
977   PetscCall(VecRestoreArray(Dgen, &dgen));
978   PetscCall(VecRestoreArray(Xgen, &xgen));
979   PetscCall(DMCompositeGather(user->dmpgrid, INSERT_VALUES, drdu_col, Dgen, Dnet));
980   PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid, &Dgen, &Dnet));
981   PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid, &Xgen, &Xnet));
982   PetscCall(VecResetArray(drdu_col));
983   PetscCall(MatDenseRestoreColumn(DRDU, &xarr));
984   PetscCall(VecDestroy(&drdu_col));
985   PetscFunctionReturn(0);
986 }
987 
988 static PetscErrorCode DRDPJacobianTranspose(TS ts, PetscReal t, Vec U, Mat drdp, Userctx *user) {
989   PetscFunctionBegin;
990   PetscFunctionReturn(0);
991 }
992 
993 PetscErrorCode ComputeSensiP(Vec lambda, Vec mu, Vec *DICDP, Userctx *user) {
994   PetscScalar *x, *y, sensip;
995   PetscInt     i;
996 
997   PetscFunctionBegin;
998   PetscCall(VecGetArray(lambda, &x));
999   PetscCall(VecGetArray(mu, &y));
1000 
1001   for (i = 0; i < 3; i++) {
1002     PetscCall(VecDot(lambda, DICDP[i], &sensip));
1003     sensip = sensip + y[i];
1004     /* PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt %" PetscInt_FMT " th parameter: %g \n",i,(double)sensip)); */
1005     y[i]   = sensip;
1006   }
1007   PetscCall(VecRestoreArray(mu, &y));
1008   PetscFunctionReturn(0);
1009 }
1010 
1011 int main(int argc, char **argv) {
1012   Userctx      user;
1013   Vec          p;
1014   PetscScalar *x_ptr;
1015   PetscMPIInt  size;
1016   PetscInt     i;
1017   PetscViewer  Xview, Ybusview;
1018   PetscInt    *idx2;
1019   Tao          tao;
1020   KSP          ksp;
1021   PC           pc;
1022   Vec          lowerb, upperb;
1023 
1024   PetscFunctionBeginUser;
1025   PetscCall(PetscInitialize(&argc, &argv, "petscoptions", help));
1026   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
1027   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs");
1028 
1029   user.jacp_flg   = PETSC_TRUE;
1030   user.neqs_gen   = 9 * ngen; /* # eqs. for generator subsystem */
1031   user.neqs_net   = 2 * nbus; /* # eqs. for network subsystem   */
1032   user.neqs_pgrid = user.neqs_gen + user.neqs_net;
1033 
1034   /* Create indices for differential and algebraic equations */
1035   PetscCall(PetscMalloc1(7 * ngen, &idx2));
1036   for (i = 0; i < ngen; i++) {
1037     idx2[7 * i]     = 9 * i;
1038     idx2[7 * i + 1] = 9 * i + 1;
1039     idx2[7 * i + 2] = 9 * i + 2;
1040     idx2[7 * i + 3] = 9 * i + 3;
1041     idx2[7 * i + 4] = 9 * i + 6;
1042     idx2[7 * i + 5] = 9 * i + 7;
1043     idx2[7 * i + 6] = 9 * i + 8;
1044   }
1045   PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, 7 * ngen, idx2, PETSC_COPY_VALUES, &user.is_diff));
1046   PetscCall(ISComplement(user.is_diff, 0, user.neqs_pgrid, &user.is_alg));
1047   PetscCall(PetscFree(idx2));
1048 
1049   /* Set run time options */
1050   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Transient stability fault options", "");
1051   {
1052     user.tfaulton  = 1.0;
1053     user.tfaultoff = 1.2;
1054     user.Rfault    = 0.0001;
1055     user.faultbus  = 8;
1056     PetscCall(PetscOptionsReal("-tfaulton", "", "", user.tfaulton, &user.tfaulton, NULL));
1057     PetscCall(PetscOptionsReal("-tfaultoff", "", "", user.tfaultoff, &user.tfaultoff, NULL));
1058     PetscCall(PetscOptionsInt("-faultbus", "", "", user.faultbus, &user.faultbus, NULL));
1059     user.t0   = 0.0;
1060     user.tmax = 1.3;
1061     PetscCall(PetscOptionsReal("-t0", "", "", user.t0, &user.t0, NULL));
1062     PetscCall(PetscOptionsReal("-tmax", "", "", user.tmax, &user.tmax, NULL));
1063     user.freq_u = 61.0;
1064     user.freq_l = 59.0;
1065     user.pow    = 2;
1066     PetscCall(PetscOptionsReal("-frequ", "", "", user.freq_u, &user.freq_u, NULL));
1067     PetscCall(PetscOptionsReal("-freql", "", "", user.freq_l, &user.freq_l, NULL));
1068     PetscCall(PetscOptionsInt("-pow", "", "", user.pow, &user.pow, NULL));
1069   }
1070   PetscOptionsEnd();
1071 
1072   /* Create DMs for generator and network subsystems */
1073   PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, user.neqs_gen, 1, 1, NULL, &user.dmgen));
1074   PetscCall(DMSetOptionsPrefix(user.dmgen, "dmgen_"));
1075   PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, user.neqs_net, 1, 1, NULL, &user.dmnet));
1076   PetscCall(DMSetOptionsPrefix(user.dmnet, "dmnet_"));
1077   PetscCall(DMSetFromOptions(user.dmnet));
1078   PetscCall(DMSetUp(user.dmnet));
1079   /* Create a composite DM packer and add the two DMs */
1080   PetscCall(DMCompositeCreate(PETSC_COMM_WORLD, &user.dmpgrid));
1081   PetscCall(DMSetOptionsPrefix(user.dmpgrid, "pgrid_"));
1082   PetscCall(DMSetFromOptions(user.dmgen));
1083   PetscCall(DMSetUp(user.dmgen));
1084   PetscCall(DMCompositeAddDM(user.dmpgrid, user.dmgen));
1085   PetscCall(DMCompositeAddDM(user.dmpgrid, user.dmnet));
1086 
1087   /* Read initial voltage vector and Ybus */
1088   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "X.bin", FILE_MODE_READ, &Xview));
1089   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "Ybus.bin", FILE_MODE_READ, &Ybusview));
1090 
1091   PetscCall(VecCreate(PETSC_COMM_WORLD, &user.V0));
1092   PetscCall(VecSetSizes(user.V0, PETSC_DECIDE, user.neqs_net));
1093   PetscCall(VecLoad(user.V0, Xview));
1094 
1095   PetscCall(MatCreate(PETSC_COMM_WORLD, &user.Ybus));
1096   PetscCall(MatSetSizes(user.Ybus, PETSC_DECIDE, PETSC_DECIDE, user.neqs_net, user.neqs_net));
1097   PetscCall(MatSetType(user.Ybus, MATBAIJ));
1098   /*  PetscCall(MatSetBlockSize(ctx->Ybus,2)); */
1099   PetscCall(MatLoad(user.Ybus, Ybusview));
1100 
1101   PetscCall(PetscViewerDestroy(&Xview));
1102   PetscCall(PetscViewerDestroy(&Ybusview));
1103 
1104   /* Allocate space for Jacobians */
1105   PetscCall(MatCreate(PETSC_COMM_WORLD, &user.J));
1106   PetscCall(MatSetSizes(user.J, PETSC_DECIDE, PETSC_DECIDE, user.neqs_pgrid, user.neqs_pgrid));
1107   PetscCall(MatSetFromOptions(user.J));
1108   PetscCall(PreallocateJacobian(user.J, &user));
1109 
1110   PetscCall(MatCreate(PETSC_COMM_WORLD, &user.Jacp));
1111   PetscCall(MatSetSizes(user.Jacp, PETSC_DECIDE, PETSC_DECIDE, user.neqs_pgrid, 3));
1112   PetscCall(MatSetFromOptions(user.Jacp));
1113   PetscCall(MatSetUp(user.Jacp));
1114   PetscCall(MatZeroEntries(user.Jacp)); /* initialize to zeros */
1115 
1116   PetscCall(MatCreateDense(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, 3, 1, NULL, &user.DRDP));
1117   PetscCall(MatSetUp(user.DRDP));
1118   PetscCall(MatCreateDense(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, user.neqs_pgrid, 1, NULL, &user.DRDU));
1119   PetscCall(MatSetUp(user.DRDU));
1120 
1121   /* Create TAO solver and set desired solution method */
1122   PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao));
1123   PetscCall(TaoSetType(tao, TAOBLMVM));
1124   /*
1125      Optimization starts
1126   */
1127   /* Set initial solution guess */
1128   PetscCall(VecCreateSeq(PETSC_COMM_WORLD, 3, &p));
1129   PetscCall(VecGetArray(p, &x_ptr));
1130   x_ptr[0] = PG[0];
1131   x_ptr[1] = PG[1];
1132   x_ptr[2] = PG[2];
1133   PetscCall(VecRestoreArray(p, &x_ptr));
1134 
1135   PetscCall(TaoSetSolution(tao, p));
1136   /* Set routine for function and gradient evaluation */
1137   PetscCall(TaoSetObjectiveAndGradient(tao, NULL, FormFunctionGradient, &user));
1138 
1139   /* Set bounds for the optimization */
1140   PetscCall(VecDuplicate(p, &lowerb));
1141   PetscCall(VecDuplicate(p, &upperb));
1142   PetscCall(VecGetArray(lowerb, &x_ptr));
1143   x_ptr[0] = 0.5;
1144   x_ptr[1] = 0.5;
1145   x_ptr[2] = 0.5;
1146   PetscCall(VecRestoreArray(lowerb, &x_ptr));
1147   PetscCall(VecGetArray(upperb, &x_ptr));
1148   x_ptr[0] = 2.0;
1149   x_ptr[1] = 2.0;
1150   x_ptr[2] = 2.0;
1151   PetscCall(VecRestoreArray(upperb, &x_ptr));
1152   PetscCall(TaoSetVariableBounds(tao, lowerb, upperb));
1153 
1154   /* Check for any TAO command line options */
1155   PetscCall(TaoSetFromOptions(tao));
1156   PetscCall(TaoGetKSP(tao, &ksp));
1157   if (ksp) {
1158     PetscCall(KSPGetPC(ksp, &pc));
1159     PetscCall(PCSetType(pc, PCNONE));
1160   }
1161 
1162   /* SOLVE THE APPLICATION */
1163   PetscCall(TaoSolve(tao));
1164 
1165   PetscCall(VecView(p, PETSC_VIEWER_STDOUT_WORLD));
1166   /* Free TAO data structures */
1167   PetscCall(TaoDestroy(&tao));
1168 
1169   PetscCall(DMDestroy(&user.dmgen));
1170   PetscCall(DMDestroy(&user.dmnet));
1171   PetscCall(DMDestroy(&user.dmpgrid));
1172   PetscCall(ISDestroy(&user.is_diff));
1173   PetscCall(ISDestroy(&user.is_alg));
1174 
1175   PetscCall(MatDestroy(&user.J));
1176   PetscCall(MatDestroy(&user.Jacp));
1177   PetscCall(MatDestroy(&user.Ybus));
1178   /* PetscCall(MatDestroy(&user.Sol)); */
1179   PetscCall(VecDestroy(&user.V0));
1180   PetscCall(VecDestroy(&p));
1181   PetscCall(VecDestroy(&lowerb));
1182   PetscCall(VecDestroy(&upperb));
1183   PetscCall(MatDestroy(&user.DRDU));
1184   PetscCall(MatDestroy(&user.DRDP));
1185   PetscCall(PetscFinalize());
1186   return 0;
1187 }
1188 
1189 /* ------------------------------------------------------------------ */
1190 /*
1191    FormFunction - Evaluates the function and corresponding gradient.
1192 
1193    Input Parameters:
1194    tao - the Tao context
1195    X   - the input vector
1196    ptr - optional user-defined context, as set by TaoSetObjectiveAndGradient()
1197 
1198    Output Parameters:
1199    f   - the newly evaluated function
1200    G   - the newly evaluated gradient
1201 */
1202 PetscErrorCode FormFunctionGradient(Tao tao, Vec P, PetscReal *f, Vec G, void *ctx0) {
1203   TS           ts, quadts;
1204   SNES         snes_alg;
1205   Userctx     *ctx = (Userctx *)ctx0;
1206   Vec          X;
1207   PetscInt     i;
1208   /* sensitivity context */
1209   PetscScalar *x_ptr;
1210   Vec          lambda[1], q;
1211   Vec          mu[1];
1212   PetscInt     steps1, steps2, steps3;
1213   Vec          DICDP[3];
1214   Vec          F_alg;
1215   PetscInt     row_loc, col_loc;
1216   PetscScalar  val;
1217   Vec          Xdot;
1218 
1219   PetscFunctionBegin;
1220   PetscCall(VecGetArrayRead(P, (const PetscScalar **)&x_ptr));
1221   PG[0] = x_ptr[0];
1222   PG[1] = x_ptr[1];
1223   PG[2] = x_ptr[2];
1224   PetscCall(VecRestoreArrayRead(P, (const PetscScalar **)&x_ptr));
1225 
1226   ctx->stepnum = 0;
1227 
1228   PetscCall(DMCreateGlobalVector(ctx->dmpgrid, &X));
1229 
1230   /* Create matrix to save solutions at each time step */
1231   /* PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,ctx->neqs_pgrid+1,1002,NULL,&ctx->Sol)); */
1232   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1233      Create timestepping solver context
1234      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1235   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
1236   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
1237   PetscCall(TSSetType(ts, TSCN));
1238   PetscCall(TSSetIFunction(ts, NULL, (TSIFunction)IFunction, ctx));
1239   PetscCall(TSSetIJacobian(ts, ctx->J, ctx->J, (TSIJacobian)IJacobian, ctx));
1240   PetscCall(TSSetApplicationContext(ts, ctx));
1241   /*   Set RHS JacobianP */
1242   PetscCall(TSSetRHSJacobianP(ts, ctx->Jacp, RHSJacobianP, ctx));
1243 
1244   PetscCall(TSCreateQuadratureTS(ts, PETSC_FALSE, &quadts));
1245   PetscCall(TSSetRHSFunction(quadts, NULL, (TSRHSFunction)CostIntegrand, ctx));
1246   PetscCall(TSSetRHSJacobian(quadts, ctx->DRDU, ctx->DRDU, (TSRHSJacobian)DRDUJacobianTranspose, ctx));
1247   PetscCall(TSSetRHSJacobianP(quadts, ctx->DRDP, (TSRHSJacobianP)DRDPJacobianTranspose, ctx));
1248 
1249   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1250      Set initial conditions
1251    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1252   PetscCall(SetInitialGuess(X, ctx));
1253 
1254   /* Approximate DICDP with finite difference, we want to zero out network variables */
1255   for (i = 0; i < 3; i++) { PetscCall(VecDuplicate(X, &DICDP[i])); }
1256   PetscCall(DICDPFiniteDifference(X, DICDP, ctx));
1257 
1258   PetscCall(VecDuplicate(X, &F_alg));
1259   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes_alg));
1260   PetscCall(SNESSetFunction(snes_alg, F_alg, AlgFunction, ctx));
1261   PetscCall(MatZeroEntries(ctx->J));
1262   PetscCall(SNESSetJacobian(snes_alg, ctx->J, ctx->J, AlgJacobian, ctx));
1263   PetscCall(SNESSetOptionsPrefix(snes_alg, "alg_"));
1264   PetscCall(SNESSetFromOptions(snes_alg));
1265   ctx->alg_flg = PETSC_TRUE;
1266   /* Solve the algebraic equations */
1267   PetscCall(SNESSolve(snes_alg, NULL, X));
1268 
1269   /* Just to set up the Jacobian structure */
1270   PetscCall(VecDuplicate(X, &Xdot));
1271   PetscCall(IJacobian(ts, 0.0, X, Xdot, 0.0, ctx->J, ctx->J, ctx));
1272   PetscCall(VecDestroy(&Xdot));
1273 
1274   ctx->stepnum++;
1275 
1276   /*
1277     Save trajectory of solution so that TSAdjointSolve() may be used
1278   */
1279   PetscCall(TSSetSaveTrajectory(ts));
1280 
1281   PetscCall(TSSetTimeStep(ts, 0.01));
1282   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
1283   PetscCall(TSSetFromOptions(ts));
1284   /* PetscCall(TSSetPostStep(ts,SaveSolution)); */
1285 
1286   /* Prefault period */
1287   ctx->alg_flg = PETSC_FALSE;
1288   PetscCall(TSSetTime(ts, 0.0));
1289   PetscCall(TSSetMaxTime(ts, ctx->tfaulton));
1290   PetscCall(TSSolve(ts, X));
1291   PetscCall(TSGetStepNumber(ts, &steps1));
1292 
1293   /* Create the nonlinear solver for solving the algebraic system */
1294   /* Note that although the algebraic system needs to be solved only for
1295      Idq and V, we reuse the entire system including xgen. The xgen
1296      variables are held constant by setting their residuals to 0 and
1297      putting a 1 on the Jacobian diagonal for xgen rows
1298   */
1299   PetscCall(MatZeroEntries(ctx->J));
1300 
1301   /* Apply disturbance - resistive fault at ctx->faultbus */
1302   /* This is done by adding shunt conductance to the diagonal location
1303      in the Ybus matrix */
1304   row_loc = 2 * ctx->faultbus;
1305   col_loc = 2 * ctx->faultbus + 1; /* Location for G */
1306   val     = 1 / ctx->Rfault;
1307   PetscCall(MatSetValues(ctx->Ybus, 1, &row_loc, 1, &col_loc, &val, ADD_VALUES));
1308   row_loc = 2 * ctx->faultbus + 1;
1309   col_loc = 2 * ctx->faultbus; /* Location for G */
1310   val     = 1 / ctx->Rfault;
1311   PetscCall(MatSetValues(ctx->Ybus, 1, &row_loc, 1, &col_loc, &val, ADD_VALUES));
1312 
1313   PetscCall(MatAssemblyBegin(ctx->Ybus, MAT_FINAL_ASSEMBLY));
1314   PetscCall(MatAssemblyEnd(ctx->Ybus, MAT_FINAL_ASSEMBLY));
1315 
1316   ctx->alg_flg = PETSC_TRUE;
1317   /* Solve the algebraic equations */
1318   PetscCall(SNESSolve(snes_alg, NULL, X));
1319 
1320   ctx->stepnum++;
1321 
1322   /* Disturbance period */
1323   ctx->alg_flg = PETSC_FALSE;
1324   PetscCall(TSSetTime(ts, ctx->tfaulton));
1325   PetscCall(TSSetMaxTime(ts, ctx->tfaultoff));
1326   PetscCall(TSSolve(ts, X));
1327   PetscCall(TSGetStepNumber(ts, &steps2));
1328   steps2 -= steps1;
1329 
1330   /* Remove the fault */
1331   row_loc = 2 * ctx->faultbus;
1332   col_loc = 2 * ctx->faultbus + 1;
1333   val     = -1 / ctx->Rfault;
1334   PetscCall(MatSetValues(ctx->Ybus, 1, &row_loc, 1, &col_loc, &val, ADD_VALUES));
1335   row_loc = 2 * ctx->faultbus + 1;
1336   col_loc = 2 * ctx->faultbus;
1337   val     = -1 / ctx->Rfault;
1338   PetscCall(MatSetValues(ctx->Ybus, 1, &row_loc, 1, &col_loc, &val, ADD_VALUES));
1339 
1340   PetscCall(MatAssemblyBegin(ctx->Ybus, MAT_FINAL_ASSEMBLY));
1341   PetscCall(MatAssemblyEnd(ctx->Ybus, MAT_FINAL_ASSEMBLY));
1342 
1343   PetscCall(MatZeroEntries(ctx->J));
1344 
1345   ctx->alg_flg = PETSC_TRUE;
1346 
1347   /* Solve the algebraic equations */
1348   PetscCall(SNESSolve(snes_alg, NULL, X));
1349 
1350   ctx->stepnum++;
1351 
1352   /* Post-disturbance period */
1353   ctx->alg_flg = PETSC_TRUE;
1354   PetscCall(TSSetTime(ts, ctx->tfaultoff));
1355   PetscCall(TSSetMaxTime(ts, ctx->tmax));
1356   PetscCall(TSSolve(ts, X));
1357   PetscCall(TSGetStepNumber(ts, &steps3));
1358   steps3 -= steps2;
1359   steps3 -= steps1;
1360 
1361   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1362      Adjoint model starts here
1363      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1364   PetscCall(TSSetPostStep(ts, NULL));
1365   PetscCall(MatCreateVecs(ctx->J, &lambda[0], NULL));
1366   /*   Set initial conditions for the adjoint integration */
1367   PetscCall(VecZeroEntries(lambda[0]));
1368 
1369   PetscCall(MatCreateVecs(ctx->Jacp, &mu[0], NULL));
1370   PetscCall(VecZeroEntries(mu[0]));
1371   PetscCall(TSSetCostGradients(ts, 1, lambda, mu));
1372 
1373   PetscCall(TSAdjointSetSteps(ts, steps3));
1374   PetscCall(TSAdjointSolve(ts));
1375 
1376   PetscCall(MatZeroEntries(ctx->J));
1377   /* Applying disturbance - resistive fault at ctx->faultbus */
1378   /* This is done by deducting shunt conductance to the diagonal location
1379      in the Ybus matrix */
1380   row_loc = 2 * ctx->faultbus;
1381   col_loc = 2 * ctx->faultbus + 1; /* Location for G */
1382   val     = 1. / ctx->Rfault;
1383   PetscCall(MatSetValues(ctx->Ybus, 1, &row_loc, 1, &col_loc, &val, ADD_VALUES));
1384   row_loc = 2 * ctx->faultbus + 1;
1385   col_loc = 2 * ctx->faultbus; /* Location for G */
1386   val     = 1. / ctx->Rfault;
1387   PetscCall(MatSetValues(ctx->Ybus, 1, &row_loc, 1, &col_loc, &val, ADD_VALUES));
1388 
1389   PetscCall(MatAssemblyBegin(ctx->Ybus, MAT_FINAL_ASSEMBLY));
1390   PetscCall(MatAssemblyEnd(ctx->Ybus, MAT_FINAL_ASSEMBLY));
1391 
1392   /*   Set number of steps for the adjoint integration */
1393   PetscCall(TSAdjointSetSteps(ts, steps2));
1394   PetscCall(TSAdjointSolve(ts));
1395 
1396   PetscCall(MatZeroEntries(ctx->J));
1397   /* remove the fault */
1398   row_loc = 2 * ctx->faultbus;
1399   col_loc = 2 * ctx->faultbus + 1; /* Location for G */
1400   val     = -1. / ctx->Rfault;
1401   PetscCall(MatSetValues(ctx->Ybus, 1, &row_loc, 1, &col_loc, &val, ADD_VALUES));
1402   row_loc = 2 * ctx->faultbus + 1;
1403   col_loc = 2 * ctx->faultbus; /* Location for G */
1404   val     = -1. / ctx->Rfault;
1405   PetscCall(MatSetValues(ctx->Ybus, 1, &row_loc, 1, &col_loc, &val, ADD_VALUES));
1406 
1407   PetscCall(MatAssemblyBegin(ctx->Ybus, MAT_FINAL_ASSEMBLY));
1408   PetscCall(MatAssemblyEnd(ctx->Ybus, MAT_FINAL_ASSEMBLY));
1409 
1410   /*   Set number of steps for the adjoint integration */
1411   PetscCall(TSAdjointSetSteps(ts, steps1));
1412   PetscCall(TSAdjointSolve(ts));
1413 
1414   PetscCall(ComputeSensiP(lambda[0], mu[0], DICDP, ctx));
1415   PetscCall(VecCopy(mu[0], G));
1416 
1417   PetscCall(TSGetQuadratureTS(ts, NULL, &quadts));
1418   PetscCall(TSGetSolution(quadts, &q));
1419   PetscCall(VecGetArray(q, &x_ptr));
1420   *f       = x_ptr[0];
1421   x_ptr[0] = 0;
1422   PetscCall(VecRestoreArray(q, &x_ptr));
1423 
1424   PetscCall(VecDestroy(&lambda[0]));
1425   PetscCall(VecDestroy(&mu[0]));
1426 
1427   PetscCall(SNESDestroy(&snes_alg));
1428   PetscCall(VecDestroy(&F_alg));
1429   PetscCall(VecDestroy(&X));
1430   PetscCall(TSDestroy(&ts));
1431   for (i = 0; i < 3; i++) { PetscCall(VecDestroy(&DICDP[i])); }
1432   PetscFunctionReturn(0);
1433 }
1434 
1435 /*TEST
1436 
1437    build:
1438       requires: double !complex !defined(PETSC_USE_64BIT_INDICES)
1439 
1440    test:
1441       args: -viewer_binary_skip_info -tao_monitor -tao_gttol .2
1442       localrunfiles: petscoptions X.bin Ybus.bin
1443 
1444 TEST*/
1445