xref: /petsc/src/ts/tutorials/power_grid/stability_9bus/ex9busoptfd.c (revision 58d68138c660dfb4e9f5b03334792cd4f2ffd7cc)
1 static char help[] = "Using finite difference for the problem in ex9busopt.c \n\n";
2 
3 /*
4   Use finite difference approximations to solve the same optimization problem as in ex9busopt.c.
5  */
6 
7 #include <petsctao.h>
8 #include <petscts.h>
9 #include <petscdm.h>
10 #include <petscdmda.h>
11 #include <petscdmcomposite.h>
12 
13 PetscErrorCode FormFunction(Tao, Vec, PetscReal *, void *);
14 
15 #define freq 60
16 #define w_s  (2 * PETSC_PI * freq)
17 
18 /* Sizes and indices */
19 const PetscInt nbus    = 9;         /* Number of network buses */
20 const PetscInt ngen    = 3;         /* Number of generators */
21 const PetscInt nload   = 3;         /* Number of loads */
22 const PetscInt gbus[3] = {0, 1, 2}; /* Buses at which generators are incident */
23 const PetscInt lbus[3] = {4, 5, 7}; /* Buses at which loads are incident */
24 
25 /* Generator real and reactive powers (found via loadflow) */
26 PetscScalar       PG[3]   = {0.69, 1.59, 0.69};
27 /* PetscScalar PG[3] = {0.716786142395021,1.630000000000000,0.850000000000000};*/
28 const PetscScalar QG[3]   = {0.270702180178785, 0.066120127797275, -0.108402221791588};
29 /* Generator constants */
30 const PetscScalar H[3]    = {23.64, 6.4, 3.01};       /* Inertia constant */
31 const PetscScalar Rs[3]   = {0.0, 0.0, 0.0};          /* Stator Resistance */
32 const PetscScalar Xd[3]   = {0.146, 0.8958, 1.3125};  /* d-axis reactance */
33 const PetscScalar Xdp[3]  = {0.0608, 0.1198, 0.1813}; /* d-axis transient reactance */
34 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 */
35 const PetscScalar Xqp[3]  = {0.0969, 0.1969, 0.25};   /* q-axis transient reactance */
36 const PetscScalar Td0p[3] = {8.96, 6.0, 5.89};        /* d-axis open circuit time constant */
37 const PetscScalar Tq0p[3] = {0.31, 0.535, 0.6};       /* q-axis open circuit time constant */
38 PetscScalar       M[3];                               /* M = 2*H/w_s */
39 PetscScalar       D[3];                               /* D = 0.1*M */
40 
41 PetscScalar       TM[3]; /* Mechanical Torque */
42 /* Exciter system constants */
43 const PetscScalar KA[3] = {20.0, 20.0, 20.0};    /* Voltage regulartor gain constant */
44 const PetscScalar TA[3] = {0.2, 0.2, 0.2};       /* Voltage regulator time constant */
45 const PetscScalar KE[3] = {1.0, 1.0, 1.0};       /* Exciter gain constant */
46 const PetscScalar TE[3] = {0.314, 0.314, 0.314}; /* Exciter time constant */
47 const PetscScalar KF[3] = {0.063, 0.063, 0.063}; /* Feedback stabilizer gain constant */
48 const PetscScalar TF[3] = {0.35, 0.35, 0.35};    /* Feedback stabilizer time constant */
49 const PetscScalar k1[3] = {0.0039, 0.0039, 0.0039};
50 const PetscScalar k2[3] = {1.555, 1.555, 1.555}; /* k1 and k2 for calculating the saturation function SE = k1*exp(k2*Efd) */
51 
52 PetscScalar       Vref[3];
53 /* Load constants
54   We use a composite load model that describes the load and reactive powers at each time instant as follows
55   P(t) = \sum\limits_{i=0}^ld_nsegsp \ld_alphap_i*P_D0(\frac{V_m(t)}{V_m0})^\ld_betap_i
56   Q(t) = \sum\limits_{i=0}^ld_nsegsq \ld_alphaq_i*Q_D0(\frac{V_m(t)}{V_m0})^\ld_betaq_i
57   where
58     ld_nsegsp,ld_nsegsq - Number of individual load models for real and reactive power loads
59     ld_alphap,ld_alphap - Percentage contribution (weights) or loads
60     P_D0                - Real power load
61     Q_D0                - Reactive power load
62     V_m(t)              - Voltage magnitude at time t
63     V_m0                - Voltage magnitude at t = 0
64     ld_betap, ld_betaq  - exponents describing the load model for real and reactive part
65 
66     Note: All loads have the same characteristic currently.
67 */
68 const PetscScalar PD0[3]       = {1.25, 0.9, 1.0};
69 const PetscScalar QD0[3]       = {0.5, 0.3, 0.35};
70 const PetscInt    ld_nsegsp[3] = {3, 3, 3};
71 const PetscScalar ld_alphap[3] = {1.0, 0.0, 0.0};
72 const PetscScalar ld_betap[3]  = {2.0, 1.0, 0.0};
73 const PetscInt    ld_nsegsq[3] = {3, 3, 3};
74 const PetscScalar ld_alphaq[3] = {1.0, 0.0, 0.0};
75 const PetscScalar ld_betaq[3]  = {2.0, 1.0, 0.0};
76 
77 typedef struct {
78   DM          dmgen, dmnet;        /* DMs to manage generator and network subsystem */
79   DM          dmpgrid;             /* Composite DM to manage the entire power grid */
80   Mat         Ybus;                /* Network admittance matrix */
81   Vec         V0;                  /* Initial voltage vector (Power flow solution) */
82   PetscReal   tfaulton, tfaultoff; /* Fault on and off times */
83   PetscInt    faultbus;            /* Fault bus */
84   PetscScalar Rfault;
85   PetscReal   t0, tmax;
86   PetscInt    neqs_gen, neqs_net, neqs_pgrid;
87   Mat         Sol; /* Matrix to save solution at each time step */
88   PetscInt    stepnum;
89   PetscBool   alg_flg;
90   PetscReal   t;
91   IS          is_diff;        /* indices for differential equations */
92   IS          is_alg;         /* indices for algebraic equations */
93   PetscReal   freq_u, freq_l; /* upper and lower frequency limit */
94   PetscInt    pow;            /* power coefficient used in the cost function */
95   Vec         vec_q;
96 } Userctx;
97 
98 /* Converts from machine frame (dq) to network (phase a real,imag) reference frame */
99 PetscErrorCode dq2ri(PetscScalar Fd, PetscScalar Fq, PetscScalar delta, PetscScalar *Fr, PetscScalar *Fi) {
100   PetscFunctionBegin;
101   *Fr = Fd * PetscSinScalar(delta) + Fq * PetscCosScalar(delta);
102   *Fi = -Fd * PetscCosScalar(delta) + Fq * PetscSinScalar(delta);
103   PetscFunctionReturn(0);
104 }
105 
106 /* Converts from network frame ([phase a real,imag) to machine (dq) reference frame */
107 PetscErrorCode ri2dq(PetscScalar Fr, PetscScalar Fi, PetscScalar delta, PetscScalar *Fd, PetscScalar *Fq) {
108   PetscFunctionBegin;
109   *Fd = Fr * PetscSinScalar(delta) - Fi * PetscCosScalar(delta);
110   *Fq = Fr * PetscCosScalar(delta) + Fi * PetscSinScalar(delta);
111   PetscFunctionReturn(0);
112 }
113 
114 PetscErrorCode SetInitialGuess(Vec X, Userctx *user) {
115   Vec          Xgen, Xnet;
116   PetscScalar *xgen, *xnet;
117   PetscInt     i, idx = 0;
118   PetscScalar  Vr, Vi, IGr, IGi, Vm, Vm2;
119   PetscScalar  Eqp, Edp, delta;
120   PetscScalar  Efd, RF, VR; /* Exciter variables */
121   PetscScalar  Id, Iq;      /* Generator dq axis currents */
122   PetscScalar  theta, Vd, Vq, SE;
123 
124   PetscFunctionBegin;
125   M[0] = 2 * H[0] / w_s;
126   M[1] = 2 * H[1] / w_s;
127   M[2] = 2 * H[2] / w_s;
128   D[0] = 0.1 * M[0];
129   D[1] = 0.1 * M[1];
130   D[2] = 0.1 * M[2];
131 
132   PetscCall(DMCompositeGetLocalVectors(user->dmpgrid, &Xgen, &Xnet));
133 
134   /* Network subsystem initialization */
135   PetscCall(VecCopy(user->V0, Xnet));
136 
137   /* Generator subsystem initialization */
138   PetscCall(VecGetArray(Xgen, &xgen));
139   PetscCall(VecGetArray(Xnet, &xnet));
140 
141   for (i = 0; i < ngen; i++) {
142     Vr  = xnet[2 * gbus[i]];     /* Real part of generator terminal voltage */
143     Vi  = xnet[2 * gbus[i] + 1]; /* Imaginary part of the generator terminal voltage */
144     Vm  = PetscSqrtScalar(Vr * Vr + Vi * Vi);
145     Vm2 = Vm * Vm;
146     IGr = (Vr * PG[i] + Vi * QG[i]) / Vm2;
147     IGi = (Vi * PG[i] - Vr * QG[i]) / Vm2;
148 
149     delta = PetscAtan2Real(Vi + Xq[i] * IGr, Vr - Xq[i] * IGi); /* Machine angle */
150 
151     theta = PETSC_PI / 2.0 - delta;
152 
153     Id = IGr * PetscCosScalar(theta) - IGi * PetscSinScalar(theta); /* d-axis stator current */
154     Iq = IGr * PetscSinScalar(theta) + IGi * PetscCosScalar(theta); /* q-axis stator current */
155 
156     Vd = Vr * PetscCosScalar(theta) - Vi * PetscSinScalar(theta);
157     Vq = Vr * PetscSinScalar(theta) + Vi * PetscCosScalar(theta);
158 
159     Edp = Vd + Rs[i] * Id - Xqp[i] * Iq; /* d-axis transient EMF */
160     Eqp = Vq + Rs[i] * Iq + Xdp[i] * Id; /* q-axis transient EMF */
161 
162     TM[i] = PG[i];
163 
164     /* The generator variables are ordered as [Eqp,Edp,delta,w,Id,Iq] */
165     xgen[idx]     = Eqp;
166     xgen[idx + 1] = Edp;
167     xgen[idx + 2] = delta;
168     xgen[idx + 3] = w_s;
169 
170     idx = idx + 4;
171 
172     xgen[idx]     = Id;
173     xgen[idx + 1] = Iq;
174 
175     idx = idx + 2;
176 
177     /* Exciter */
178     Efd = Eqp + (Xd[i] - Xdp[i]) * Id;
179     SE  = k1[i] * PetscExpScalar(k2[i] * Efd);
180     VR  = KE[i] * Efd + SE;
181     RF  = KF[i] * Efd / TF[i];
182 
183     xgen[idx]     = Efd;
184     xgen[idx + 1] = RF;
185     xgen[idx + 2] = VR;
186 
187     Vref[i] = Vm + (VR / KA[i]);
188 
189     idx = idx + 3;
190   }
191 
192   PetscCall(VecRestoreArray(Xgen, &xgen));
193   PetscCall(VecRestoreArray(Xnet, &xnet));
194 
195   /* PetscCall(VecView(Xgen,0)); */
196   PetscCall(DMCompositeGather(user->dmpgrid, INSERT_VALUES, X, Xgen, Xnet));
197   PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid, &Xgen, &Xnet));
198   PetscFunctionReturn(0);
199 }
200 
201 /* Computes F = [-f(x,y);g(x,y)] */
202 PetscErrorCode ResidualFunction(SNES snes, Vec X, Vec F, Userctx *user) {
203   Vec          Xgen, Xnet, Fgen, Fnet;
204   PetscScalar *xgen, *xnet, *fgen, *fnet;
205   PetscInt     i, idx = 0;
206   PetscScalar  Vr, Vi, Vm, Vm2;
207   PetscScalar  Eqp, Edp, delta, w; /* Generator variables */
208   PetscScalar  Efd, RF, VR;        /* Exciter variables */
209   PetscScalar  Id, Iq;             /* Generator dq axis currents */
210   PetscScalar  Vd, Vq, SE;
211   PetscScalar  IGr, IGi, IDr, IDi;
212   PetscScalar  Zdq_inv[4], det;
213   PetscScalar  PD, QD, Vm0, *v0;
214   PetscInt     k;
215 
216   PetscFunctionBegin;
217   PetscCall(VecZeroEntries(F));
218   PetscCall(DMCompositeGetLocalVectors(user->dmpgrid, &Xgen, &Xnet));
219   PetscCall(DMCompositeGetLocalVectors(user->dmpgrid, &Fgen, &Fnet));
220   PetscCall(DMCompositeScatter(user->dmpgrid, X, Xgen, Xnet));
221   PetscCall(DMCompositeScatter(user->dmpgrid, F, Fgen, Fnet));
222 
223   /* Network current balance residual IG + Y*V + IL = 0. Only YV is added here.
224      The generator current injection, IG, and load current injection, ID are added later
225   */
226   /* Note that the values in Ybus are stored assuming the imaginary current balance
227      equation is ordered first followed by real current balance equation for each bus.
228      Thus imaginary current contribution goes in location 2*i, and
229      real current contribution in 2*i+1
230   */
231   PetscCall(MatMult(user->Ybus, Xnet, Fnet));
232 
233   PetscCall(VecGetArray(Xgen, &xgen));
234   PetscCall(VecGetArray(Xnet, &xnet));
235   PetscCall(VecGetArray(Fgen, &fgen));
236   PetscCall(VecGetArray(Fnet, &fnet));
237 
238   /* Generator subsystem */
239   for (i = 0; i < ngen; i++) {
240     Eqp   = xgen[idx];
241     Edp   = xgen[idx + 1];
242     delta = xgen[idx + 2];
243     w     = xgen[idx + 3];
244     Id    = xgen[idx + 4];
245     Iq    = xgen[idx + 5];
246     Efd   = xgen[idx + 6];
247     RF    = xgen[idx + 7];
248     VR    = xgen[idx + 8];
249 
250     /* Generator differential equations */
251     fgen[idx]     = (Eqp + (Xd[i] - Xdp[i]) * Id - Efd) / Td0p[i];
252     fgen[idx + 1] = (Edp - (Xq[i] - Xqp[i]) * Iq) / Tq0p[i];
253     fgen[idx + 2] = -w + w_s;
254     fgen[idx + 3] = (-TM[i] + Edp * Id + Eqp * Iq + (Xqp[i] - Xdp[i]) * Id * Iq + D[i] * (w - w_s)) / M[i];
255 
256     Vr = xnet[2 * gbus[i]];     /* Real part of generator terminal voltage */
257     Vi = xnet[2 * gbus[i] + 1]; /* Imaginary part of the generator terminal voltage */
258 
259     PetscCall(ri2dq(Vr, Vi, delta, &Vd, &Vq));
260     /* Algebraic equations for stator currents */
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 
588     dSE_dEfd = k1[i] * k2[i] * PetscExpScalar(k2[i] * Efd);
589 
590     row[0] = idx + 6;
591     col[0] = idx + 6;
592     col[1] = idx + 8;
593     val[0] = (KE[i] + dSE_dEfd) / TE[i];
594     val[1] = -1 / TE[i];
595     PetscCall(MatSetValues(J, 1, row, 2, col, val, INSERT_VALUES));
596 
597     /* Exciter differential equations */
598 
599     /*    fgen[idx+7] = (RF - KF[i]*Efd/TF[i])/TF[i]; */
600     row[0] = idx + 7;
601     col[0] = idx + 6;
602     col[1] = idx + 7;
603     val[0] = (-KF[i] / TF[i]) / TF[i];
604     val[1] = 1 / TF[i];
605     PetscCall(MatSetValues(J, 1, row, 2, col, val, INSERT_VALUES));
606 
607     /*    fgen[idx+8] = (VR - KA[i]*RF + KA[i]*KF[i]*Efd/TF[i] - KA[i]*(Vref[i] - Vm))/TA[i]; */
608     /* Vm = (Vd^2 + Vq^2)^0.5; */
609     dVm_dVd = Vd / Vm;
610     dVm_dVq = Vq / Vm;
611     dVm_dVr = dVm_dVd * dVd_dVr + dVm_dVq * dVq_dVr;
612     dVm_dVi = dVm_dVd * dVd_dVi + dVm_dVq * dVq_dVi;
613     row[0]  = idx + 8;
614     col[0]  = idx + 6;
615     col[1]  = idx + 7;
616     col[2]  = idx + 8;
617     val[0]  = (KA[i] * KF[i] / TF[i]) / TA[i];
618     val[1]  = -KA[i] / TA[i];
619     val[2]  = 1 / TA[i];
620     col[3]  = net_start + 2 * gbus[i];
621     col[4]  = net_start + 2 * gbus[i] + 1;
622     val[3]  = KA[i] * dVm_dVr / TA[i];
623     val[4]  = KA[i] * dVm_dVi / TA[i];
624     PetscCall(MatSetValues(J, 1, row, 5, col, val, INSERT_VALUES));
625     idx = idx + 9;
626   }
627 
628   for (i = 0; i < nbus; i++) {
629     PetscCall(MatGetRow(user->Ybus, 2 * i, &ncols, &cols, &yvals));
630     row[0] = net_start + 2 * i;
631     for (k = 0; k < ncols; k++) {
632       col[k] = net_start + cols[k];
633       val[k] = yvals[k];
634     }
635     PetscCall(MatSetValues(J, 1, row, ncols, col, val, INSERT_VALUES));
636     PetscCall(MatRestoreRow(user->Ybus, 2 * i, &ncols, &cols, &yvals));
637 
638     PetscCall(MatGetRow(user->Ybus, 2 * i + 1, &ncols, &cols, &yvals));
639     row[0] = net_start + 2 * i + 1;
640     for (k = 0; k < ncols; k++) {
641       col[k] = net_start + cols[k];
642       val[k] = yvals[k];
643     }
644     PetscCall(MatSetValues(J, 1, row, ncols, col, val, INSERT_VALUES));
645     PetscCall(MatRestoreRow(user->Ybus, 2 * i + 1, &ncols, &cols, &yvals));
646   }
647 
648   PetscCall(MatAssemblyBegin(J, MAT_FLUSH_ASSEMBLY));
649   PetscCall(MatAssemblyEnd(J, MAT_FLUSH_ASSEMBLY));
650 
651   PetscCall(VecGetArray(user->V0, &v0));
652   for (i = 0; i < nload; i++) {
653     Vr  = xnet[2 * lbus[i]];     /* Real part of load bus voltage */
654     Vi  = xnet[2 * lbus[i] + 1]; /* Imaginary part of the load bus voltage */
655     Vm  = PetscSqrtScalar(Vr * Vr + Vi * Vi);
656     Vm2 = Vm * Vm;
657     Vm4 = Vm2 * Vm2;
658     Vm0 = PetscSqrtScalar(v0[2 * lbus[i]] * v0[2 * lbus[i]] + v0[2 * lbus[i] + 1] * v0[2 * lbus[i] + 1]);
659     PD = QD = 0.0;
660     dPD_dVr = dPD_dVi = dQD_dVr = dQD_dVi = 0.0;
661     for (k = 0; k < ld_nsegsp[i]; k++) {
662       PD += ld_alphap[k] * PD0[i] * PetscPowScalar((Vm / Vm0), ld_betap[k]);
663       dPD_dVr += ld_alphap[k] * ld_betap[k] * PD0[i] * PetscPowScalar((1 / Vm0), ld_betap[k]) * Vr * PetscPowScalar(Vm, (ld_betap[k] - 2));
664       dPD_dVi += ld_alphap[k] * ld_betap[k] * PD0[i] * PetscPowScalar((1 / Vm0), ld_betap[k]) * Vi * PetscPowScalar(Vm, (ld_betap[k] - 2));
665     }
666     for (k = 0; k < ld_nsegsq[i]; k++) {
667       QD += ld_alphaq[k] * QD0[i] * PetscPowScalar((Vm / Vm0), ld_betaq[k]);
668       dQD_dVr += ld_alphaq[k] * ld_betaq[k] * QD0[i] * PetscPowScalar((1 / Vm0), ld_betaq[k]) * Vr * PetscPowScalar(Vm, (ld_betaq[k] - 2));
669       dQD_dVi += ld_alphaq[k] * ld_betaq[k] * QD0[i] * PetscPowScalar((1 / Vm0), ld_betaq[k]) * Vi * PetscPowScalar(Vm, (ld_betaq[k] - 2));
670     }
671 
672     /*    IDr = (PD*Vr + QD*Vi)/Vm2; */
673     /*    IDi = (-QD*Vr + PD*Vi)/Vm2; */
674 
675     dIDr_dVr = (dPD_dVr * Vr + dQD_dVr * Vi + PD) / Vm2 - ((PD * Vr + QD * Vi) * 2 * Vr) / Vm4;
676     dIDr_dVi = (dPD_dVi * Vr + dQD_dVi * Vi + QD) / Vm2 - ((PD * Vr + QD * Vi) * 2 * Vi) / Vm4;
677 
678     dIDi_dVr = (-dQD_dVr * Vr + dPD_dVr * Vi - QD) / Vm2 - ((-QD * Vr + PD * Vi) * 2 * Vr) / Vm4;
679     dIDi_dVi = (-dQD_dVi * Vr + dPD_dVi * Vi + PD) / Vm2 - ((-QD * Vr + PD * Vi) * 2 * Vi) / Vm4;
680 
681     /*    fnet[2*lbus[i]]   += IDi; */
682     row[0] = net_start + 2 * lbus[i];
683     col[0] = net_start + 2 * lbus[i];
684     col[1] = net_start + 2 * lbus[i] + 1;
685     val[0] = dIDi_dVr;
686     val[1] = dIDi_dVi;
687     PetscCall(MatSetValues(J, 1, row, 2, col, val, ADD_VALUES));
688     /*    fnet[2*lbus[i]+1] += IDr; */
689     row[0] = net_start + 2 * lbus[i] + 1;
690     col[0] = net_start + 2 * lbus[i];
691     col[1] = net_start + 2 * lbus[i] + 1;
692     val[0] = dIDr_dVr;
693     val[1] = dIDr_dVi;
694     PetscCall(MatSetValues(J, 1, row, 2, col, val, ADD_VALUES));
695   }
696   PetscCall(VecRestoreArray(user->V0, &v0));
697 
698   PetscCall(VecRestoreArray(Xgen, &xgen));
699   PetscCall(VecRestoreArray(Xnet, &xnet));
700 
701   PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid, &Xgen, &Xnet));
702 
703   PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
704   PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
705   PetscFunctionReturn(0);
706 }
707 
708 /*
709    J = [I, 0
710         dg_dx, dg_dy]
711 */
712 PetscErrorCode AlgJacobian(SNES snes, Vec X, Mat A, Mat B, void *ctx) {
713   Userctx *user = (Userctx *)ctx;
714 
715   PetscFunctionBegin;
716   PetscCall(ResidualJacobian(snes, X, A, B, ctx));
717   PetscCall(MatSetOption(A, MAT_KEEP_NONZERO_PATTERN, PETSC_TRUE));
718   PetscCall(MatZeroRowsIS(A, user->is_diff, 1.0, NULL, NULL));
719   PetscFunctionReturn(0);
720 }
721 
722 /*
723    J = [a*I-df_dx, -df_dy
724         dg_dx, dg_dy]
725 */
726 
727 PetscErrorCode IJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat A, Mat B, Userctx *user) {
728   SNES        snes;
729   PetscScalar atmp = (PetscScalar)a;
730   PetscInt    i, row;
731 
732   PetscFunctionBegin;
733   user->t = t;
734 
735   PetscCall(TSGetSNES(ts, &snes));
736   PetscCall(ResidualJacobian(snes, X, A, B, user));
737   for (i = 0; i < ngen; i++) {
738     row = 9 * i;
739     PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES));
740     row = 9 * i + 1;
741     PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES));
742     row = 9 * i + 2;
743     PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES));
744     row = 9 * i + 3;
745     PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES));
746     row = 9 * i + 6;
747     PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES));
748     row = 9 * i + 7;
749     PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES));
750     row = 9 * i + 8;
751     PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES));
752   }
753   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
754   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
755   PetscFunctionReturn(0);
756 }
757 
758 static PetscErrorCode CostIntegrand(TS ts, PetscReal t, Vec U, Vec R, Userctx *user) {
759   PetscScalar       *r;
760   const PetscScalar *u;
761   PetscInt           idx;
762   Vec                Xgen, Xnet;
763   PetscScalar       *xgen;
764   PetscInt           i;
765 
766   PetscFunctionBegin;
767   PetscCall(DMCompositeGetLocalVectors(user->dmpgrid, &Xgen, &Xnet));
768   PetscCall(DMCompositeScatter(user->dmpgrid, U, Xgen, Xnet));
769 
770   PetscCall(VecGetArray(Xgen, &xgen));
771 
772   PetscCall(VecGetArrayRead(U, &u));
773   PetscCall(VecGetArray(R, &r));
774   r[0] = 0.;
775 
776   idx = 0;
777   for (i = 0; i < ngen; i++) {
778     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);
779     idx += 9;
780   }
781   PetscCall(VecRestoreArray(R, &r));
782   PetscCall(VecRestoreArrayRead(U, &u));
783   PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid, &Xgen, &Xnet));
784   PetscFunctionReturn(0);
785 }
786 
787 static PetscErrorCode MonitorUpdateQ(TS ts, PetscInt stepnum, PetscReal time, Vec X, void *ctx0) {
788   Vec       C, *Y;
789   PetscInt  Nr;
790   PetscReal h, theta;
791   Userctx  *ctx = (Userctx *)ctx0;
792 
793   PetscFunctionBegin;
794   theta = 0.5;
795   PetscCall(TSGetStages(ts, &Nr, &Y));
796   PetscCall(TSGetTimeStep(ts, &h));
797   PetscCall(VecDuplicate(ctx->vec_q, &C));
798   /* compute integrals */
799   if (stepnum > 0) {
800     PetscCall(CostIntegrand(ts, time, X, C, ctx));
801     PetscCall(VecAXPY(ctx->vec_q, h * theta, C));
802     PetscCall(CostIntegrand(ts, time + h * theta, Y[0], C, ctx));
803     PetscCall(VecAXPY(ctx->vec_q, h * (1 - theta), C));
804   }
805   PetscCall(VecDestroy(&C));
806   PetscFunctionReturn(0);
807 }
808 
809 int main(int argc, char **argv) {
810   Userctx      user;
811   Vec          p;
812   PetscScalar *x_ptr;
813   PetscMPIInt  size;
814   PetscInt     i;
815   KSP          ksp;
816   PC           pc;
817   PetscInt    *idx2;
818   Tao          tao;
819   Vec          lowerb, upperb;
820 
821   PetscFunctionBeginUser;
822   PetscFunctionBeginUser;
823   PetscCall(PetscInitialize(&argc, &argv, "petscoptions", help));
824   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
825   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs");
826 
827   PetscCall(VecCreateSeq(PETSC_COMM_WORLD, 1, &user.vec_q));
828 
829   user.neqs_gen   = 9 * ngen; /* # eqs. for generator subsystem */
830   user.neqs_net   = 2 * nbus; /* # eqs. for network subsystem   */
831   user.neqs_pgrid = user.neqs_gen + user.neqs_net;
832 
833   /* Create indices for differential and algebraic equations */
834   PetscCall(PetscMalloc1(7 * ngen, &idx2));
835   for (i = 0; i < ngen; i++) {
836     idx2[7 * i]     = 9 * i;
837     idx2[7 * i + 1] = 9 * i + 1;
838     idx2[7 * i + 2] = 9 * i + 2;
839     idx2[7 * i + 3] = 9 * i + 3;
840     idx2[7 * i + 4] = 9 * i + 6;
841     idx2[7 * i + 5] = 9 * i + 7;
842     idx2[7 * i + 6] = 9 * i + 8;
843   }
844   PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, 7 * ngen, idx2, PETSC_COPY_VALUES, &user.is_diff));
845   PetscCall(ISComplement(user.is_diff, 0, user.neqs_pgrid, &user.is_alg));
846   PetscCall(PetscFree(idx2));
847 
848   /* Set run time options */
849   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Transient stability fault options", "");
850   {
851     user.tfaulton  = 1.0;
852     user.tfaultoff = 1.2;
853     user.Rfault    = 0.0001;
854     user.faultbus  = 8;
855     PetscCall(PetscOptionsReal("-tfaulton", "", "", user.tfaulton, &user.tfaulton, NULL));
856     PetscCall(PetscOptionsReal("-tfaultoff", "", "", user.tfaultoff, &user.tfaultoff, NULL));
857     PetscCall(PetscOptionsInt("-faultbus", "", "", user.faultbus, &user.faultbus, NULL));
858     user.t0   = 0.0;
859     user.tmax = 1.5;
860     PetscCall(PetscOptionsReal("-t0", "", "", user.t0, &user.t0, NULL));
861     PetscCall(PetscOptionsReal("-tmax", "", "", user.tmax, &user.tmax, NULL));
862     user.freq_u = 61.0;
863     user.freq_l = 59.0;
864     user.pow    = 2;
865     PetscCall(PetscOptionsReal("-frequ", "", "", user.freq_u, &user.freq_u, NULL));
866     PetscCall(PetscOptionsReal("-freql", "", "", user.freq_l, &user.freq_l, NULL));
867     PetscCall(PetscOptionsInt("-pow", "", "", user.pow, &user.pow, NULL));
868   }
869   PetscOptionsEnd();
870 
871   /* Create DMs for generator and network subsystems */
872   PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, user.neqs_gen, 1, 1, NULL, &user.dmgen));
873   PetscCall(DMSetOptionsPrefix(user.dmgen, "dmgen_"));
874   PetscCall(DMSetFromOptions(user.dmgen));
875   PetscCall(DMSetUp(user.dmgen));
876   PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, user.neqs_net, 1, 1, NULL, &user.dmnet));
877   PetscCall(DMSetOptionsPrefix(user.dmnet, "dmnet_"));
878   PetscCall(DMSetFromOptions(user.dmnet));
879   PetscCall(DMSetUp(user.dmnet));
880   /* Create a composite DM packer and add the two DMs */
881   PetscCall(DMCompositeCreate(PETSC_COMM_WORLD, &user.dmpgrid));
882   PetscCall(DMSetOptionsPrefix(user.dmpgrid, "pgrid_"));
883   PetscCall(DMCompositeAddDM(user.dmpgrid, user.dmgen));
884   PetscCall(DMCompositeAddDM(user.dmpgrid, user.dmnet));
885 
886   /* Create TAO solver and set desired solution method */
887   PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao));
888   PetscCall(TaoSetType(tao, TAOBLMVM));
889   /*
890      Optimization starts
891   */
892   /* Set initial solution guess */
893   PetscCall(VecCreateSeq(PETSC_COMM_WORLD, 3, &p));
894   PetscCall(VecGetArray(p, &x_ptr));
895   x_ptr[0] = PG[0];
896   x_ptr[1] = PG[1];
897   x_ptr[2] = PG[2];
898   PetscCall(VecRestoreArray(p, &x_ptr));
899 
900   PetscCall(TaoSetSolution(tao, p));
901   /* Set routine for function and gradient evaluation */
902   PetscCall(TaoSetObjective(tao, FormFunction, (void *)&user));
903   PetscCall(TaoSetGradient(tao, NULL, TaoDefaultComputeGradient, (void *)&user));
904 
905   /* Set bounds for the optimization */
906   PetscCall(VecDuplicate(p, &lowerb));
907   PetscCall(VecDuplicate(p, &upperb));
908   PetscCall(VecGetArray(lowerb, &x_ptr));
909   x_ptr[0] = 0.5;
910   x_ptr[1] = 0.5;
911   x_ptr[2] = 0.5;
912   PetscCall(VecRestoreArray(lowerb, &x_ptr));
913   PetscCall(VecGetArray(upperb, &x_ptr));
914   x_ptr[0] = 2.0;
915   x_ptr[1] = 2.0;
916   x_ptr[2] = 2.0;
917   PetscCall(VecRestoreArray(upperb, &x_ptr));
918   PetscCall(TaoSetVariableBounds(tao, lowerb, upperb));
919 
920   /* Check for any TAO command line options */
921   PetscCall(TaoSetFromOptions(tao));
922   PetscCall(TaoGetKSP(tao, &ksp));
923   if (ksp) {
924     PetscCall(KSPGetPC(ksp, &pc));
925     PetscCall(PCSetType(pc, PCNONE));
926   }
927 
928   /* SOLVE THE APPLICATION */
929   PetscCall(TaoSolve(tao));
930 
931   PetscCall(VecView(p, PETSC_VIEWER_STDOUT_WORLD));
932   /* Free TAO data structures */
933   PetscCall(TaoDestroy(&tao));
934   PetscCall(VecDestroy(&user.vec_q));
935   PetscCall(VecDestroy(&lowerb));
936   PetscCall(VecDestroy(&upperb));
937   PetscCall(VecDestroy(&p));
938   PetscCall(DMDestroy(&user.dmgen));
939   PetscCall(DMDestroy(&user.dmnet));
940   PetscCall(DMDestroy(&user.dmpgrid));
941   PetscCall(ISDestroy(&user.is_diff));
942   PetscCall(ISDestroy(&user.is_alg));
943   PetscCall(PetscFinalize());
944   return 0;
945 }
946 
947 /* ------------------------------------------------------------------ */
948 /*
949    FormFunction - Evaluates the function and corresponding gradient.
950 
951    Input Parameters:
952    tao - the Tao context
953    X   - the input vector
954    ptr - optional user-defined context, as set by TaoSetObjectiveAndGradient()
955 
956    Output Parameters:
957    f   - the newly evaluated function
958 */
959 PetscErrorCode FormFunction(Tao tao, Vec P, PetscReal *f, void *ctx0) {
960   TS           ts;
961   SNES         snes_alg;
962   Userctx     *ctx = (Userctx *)ctx0;
963   Vec          X;
964   Mat          J;
965   /* sensitivity context */
966   PetscScalar *x_ptr;
967   PetscViewer  Xview, Ybusview;
968   Vec          F_alg;
969   Vec          Xdot;
970   PetscInt     row_loc, col_loc;
971   PetscScalar  val;
972 
973   PetscCall(VecGetArrayRead(P, (const PetscScalar **)&x_ptr));
974   PG[0] = x_ptr[0];
975   PG[1] = x_ptr[1];
976   PG[2] = x_ptr[2];
977   PetscCall(VecRestoreArrayRead(P, (const PetscScalar **)&x_ptr));
978 
979   ctx->stepnum = 0;
980 
981   PetscCall(VecZeroEntries(ctx->vec_q));
982 
983   /* Read initial voltage vector and Ybus */
984   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "X.bin", FILE_MODE_READ, &Xview));
985   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "Ybus.bin", FILE_MODE_READ, &Ybusview));
986 
987   PetscCall(VecCreate(PETSC_COMM_WORLD, &ctx->V0));
988   PetscCall(VecSetSizes(ctx->V0, PETSC_DECIDE, ctx->neqs_net));
989   PetscCall(VecLoad(ctx->V0, Xview));
990 
991   PetscCall(MatCreate(PETSC_COMM_WORLD, &ctx->Ybus));
992   PetscCall(MatSetSizes(ctx->Ybus, PETSC_DECIDE, PETSC_DECIDE, ctx->neqs_net, ctx->neqs_net));
993   PetscCall(MatSetType(ctx->Ybus, MATBAIJ));
994   /*  PetscCall(MatSetBlockSize(ctx->Ybus,2)); */
995   PetscCall(MatLoad(ctx->Ybus, Ybusview));
996 
997   PetscCall(PetscViewerDestroy(&Xview));
998   PetscCall(PetscViewerDestroy(&Ybusview));
999 
1000   PetscCall(DMCreateGlobalVector(ctx->dmpgrid, &X));
1001 
1002   PetscCall(MatCreate(PETSC_COMM_WORLD, &J));
1003   PetscCall(MatSetSizes(J, PETSC_DECIDE, PETSC_DECIDE, ctx->neqs_pgrid, ctx->neqs_pgrid));
1004   PetscCall(MatSetFromOptions(J));
1005   PetscCall(PreallocateJacobian(J, ctx));
1006 
1007   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1008      Create timestepping solver context
1009      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1010   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
1011   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
1012   PetscCall(TSSetType(ts, TSCN));
1013   PetscCall(TSSetIFunction(ts, NULL, (TSIFunction)IFunction, ctx));
1014   PetscCall(TSSetIJacobian(ts, J, J, (TSIJacobian)IJacobian, ctx));
1015   PetscCall(TSSetApplicationContext(ts, ctx));
1016 
1017   PetscCall(TSMonitorSet(ts, MonitorUpdateQ, ctx, NULL));
1018   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1019      Set initial conditions
1020    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1021   PetscCall(SetInitialGuess(X, ctx));
1022 
1023   PetscCall(VecDuplicate(X, &F_alg));
1024   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes_alg));
1025   PetscCall(SNESSetFunction(snes_alg, F_alg, AlgFunction, ctx));
1026   PetscCall(MatZeroEntries(J));
1027   PetscCall(SNESSetJacobian(snes_alg, J, J, AlgJacobian, ctx));
1028   PetscCall(SNESSetOptionsPrefix(snes_alg, "alg_"));
1029   PetscCall(SNESSetFromOptions(snes_alg));
1030   ctx->alg_flg = PETSC_TRUE;
1031   /* Solve the algebraic equations */
1032   PetscCall(SNESSolve(snes_alg, NULL, X));
1033 
1034   /* Just to set up the Jacobian structure */
1035   PetscCall(VecDuplicate(X, &Xdot));
1036   PetscCall(IJacobian(ts, 0.0, X, Xdot, 0.0, J, J, ctx));
1037   PetscCall(VecDestroy(&Xdot));
1038 
1039   ctx->stepnum++;
1040 
1041   PetscCall(TSSetTimeStep(ts, 0.01));
1042   PetscCall(TSSetMaxTime(ts, ctx->tmax));
1043   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
1044   PetscCall(TSSetFromOptions(ts));
1045 
1046   /* Prefault period */
1047   ctx->alg_flg = PETSC_FALSE;
1048   PetscCall(TSSetTime(ts, 0.0));
1049   PetscCall(TSSetMaxTime(ts, ctx->tfaulton));
1050   PetscCall(TSSolve(ts, X));
1051 
1052   /* Create the nonlinear solver for solving the algebraic system */
1053   /* Note that although the algebraic system needs to be solved only for
1054      Idq and V, we reuse the entire system including xgen. The xgen
1055      variables are held constant by setting their residuals to 0 and
1056      putting a 1 on the Jacobian diagonal for xgen rows
1057   */
1058   PetscCall(MatZeroEntries(J));
1059 
1060   /* Apply disturbance - resistive fault at ctx->faultbus */
1061   /* This is done by adding shunt conductance to the diagonal location
1062      in the Ybus matrix */
1063   row_loc = 2 * ctx->faultbus;
1064   col_loc = 2 * ctx->faultbus + 1; /* Location for G */
1065   val     = 1 / ctx->Rfault;
1066   PetscCall(MatSetValues(ctx->Ybus, 1, &row_loc, 1, &col_loc, &val, ADD_VALUES));
1067   row_loc = 2 * ctx->faultbus + 1;
1068   col_loc = 2 * ctx->faultbus; /* Location for G */
1069   val     = 1 / ctx->Rfault;
1070   PetscCall(MatSetValues(ctx->Ybus, 1, &row_loc, 1, &col_loc, &val, ADD_VALUES));
1071 
1072   PetscCall(MatAssemblyBegin(ctx->Ybus, MAT_FINAL_ASSEMBLY));
1073   PetscCall(MatAssemblyEnd(ctx->Ybus, MAT_FINAL_ASSEMBLY));
1074 
1075   ctx->alg_flg = PETSC_TRUE;
1076   /* Solve the algebraic equations */
1077   PetscCall(SNESSolve(snes_alg, NULL, X));
1078 
1079   ctx->stepnum++;
1080 
1081   /* Disturbance period */
1082   ctx->alg_flg = PETSC_FALSE;
1083   PetscCall(TSSetTime(ts, ctx->tfaulton));
1084   PetscCall(TSSetMaxTime(ts, ctx->tfaultoff));
1085   PetscCall(TSSolve(ts, X));
1086 
1087   /* Remove the fault */
1088   row_loc = 2 * ctx->faultbus;
1089   col_loc = 2 * ctx->faultbus + 1;
1090   val     = -1 / ctx->Rfault;
1091   PetscCall(MatSetValues(ctx->Ybus, 1, &row_loc, 1, &col_loc, &val, ADD_VALUES));
1092   row_loc = 2 * ctx->faultbus + 1;
1093   col_loc = 2 * ctx->faultbus;
1094   val     = -1 / ctx->Rfault;
1095   PetscCall(MatSetValues(ctx->Ybus, 1, &row_loc, 1, &col_loc, &val, ADD_VALUES));
1096 
1097   PetscCall(MatAssemblyBegin(ctx->Ybus, MAT_FINAL_ASSEMBLY));
1098   PetscCall(MatAssemblyEnd(ctx->Ybus, MAT_FINAL_ASSEMBLY));
1099 
1100   PetscCall(MatZeroEntries(J));
1101 
1102   ctx->alg_flg = PETSC_TRUE;
1103 
1104   /* Solve the algebraic equations */
1105   PetscCall(SNESSolve(snes_alg, NULL, X));
1106 
1107   ctx->stepnum++;
1108 
1109   /* Post-disturbance period */
1110   ctx->alg_flg = PETSC_TRUE;
1111   PetscCall(TSSetTime(ts, ctx->tfaultoff));
1112   PetscCall(TSSetMaxTime(ts, ctx->tmax));
1113   PetscCall(TSSolve(ts, X));
1114 
1115   PetscCall(VecGetArray(ctx->vec_q, &x_ptr));
1116   *f = x_ptr[0];
1117   PetscCall(VecRestoreArray(ctx->vec_q, &x_ptr));
1118 
1119   PetscCall(MatDestroy(&ctx->Ybus));
1120   PetscCall(VecDestroy(&ctx->V0));
1121   PetscCall(SNESDestroy(&snes_alg));
1122   PetscCall(VecDestroy(&F_alg));
1123   PetscCall(MatDestroy(&J));
1124   PetscCall(VecDestroy(&X));
1125   PetscCall(TSDestroy(&ts));
1126 
1127   return 0;
1128 }
1129 
1130 /*TEST
1131 
1132   build:
1133       requires: double !complex !defined(USE_64BIT_INDICES)
1134 
1135    test:
1136       args: -viewer_binary_skip_info -tao_monitor -tao_gttol .2
1137       localrunfiles: petscoptions X.bin Ybus.bin
1138 
1139 TEST*/
1140