xref: /petsc/src/ts/tutorials/power_grid/stability_9bus/ex9busopt.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
1c4762a1bSJed Brown static char help[] = "Application of adjoint sensitivity analysis for power grid stability analysis of WECC 9 bus system.\n\
2c4762a1bSJed Brown This example is based on the 9-bus (node) example given in the book Power\n\
3c4762a1bSJed Brown Systems Dynamics and Stability (Chapter 7) by P. Sauer and M. A. Pai.\n\
4c4762a1bSJed Brown The power grid in this example consists of 9 buses (nodes), 3 generators,\n\
5c4762a1bSJed Brown 3 loads, and 9 transmission lines. The network equations are written\n\
6a5b23f4aSJose E. Roman in current balance form using rectangular coordinates.\n\n";
7c4762a1bSJed Brown 
8c4762a1bSJed Brown /*
9c4762a1bSJed Brown   This code demonstrates how to solve a DAE-constrained optimization problem with TAO, TSAdjoint and TS.
10c4762a1bSJed Brown   The objectivie is to find optimal parameter PG for each generator to minizie the frequency violations due to faults.
11c4762a1bSJed Brown   The problem features discontinuities and a cost function in integral form.
12c4762a1bSJed Brown   The gradient is computed with the discrete adjoint of an implicit theta method, see ex9busadj.c for details.
13c4762a1bSJed Brown */
14c4762a1bSJed Brown 
15c4762a1bSJed Brown #include <petsctao.h>
16c4762a1bSJed Brown #include <petscts.h>
17c4762a1bSJed Brown #include <petscdm.h>
18c4762a1bSJed Brown #include <petscdmda.h>
19c4762a1bSJed Brown #include <petscdmcomposite.h>
20c4762a1bSJed Brown #include <petsctime.h>
21c4762a1bSJed Brown 
22c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao, Vec, PetscReal *, Vec, void *);
23c4762a1bSJed Brown 
24c4762a1bSJed Brown #define freq 60
25c4762a1bSJed Brown #define w_s  (2 * PETSC_PI * freq)
26c4762a1bSJed Brown 
27c4762a1bSJed Brown /* Sizes and indices */
28c4762a1bSJed Brown const PetscInt nbus    = 9;         /* Number of network buses */
29c4762a1bSJed Brown const PetscInt ngen    = 3;         /* Number of generators */
30c4762a1bSJed Brown const PetscInt nload   = 3;         /* Number of loads */
31c4762a1bSJed Brown const PetscInt gbus[3] = {0, 1, 2}; /* Buses at which generators are incident */
32c4762a1bSJed Brown const PetscInt lbus[3] = {4, 5, 7}; /* Buses at which loads are incident */
33c4762a1bSJed Brown 
34c4762a1bSJed Brown /* Generator real and reactive powers (found via loadflow) */
35c4762a1bSJed Brown PetscScalar PG[3] = {0.69, 1.59, 0.69};
36c4762a1bSJed Brown /* PetscScalar PG[3] = {0.716786142395021,1.630000000000000,0.850000000000000};*/
37c4762a1bSJed Brown 
38c4762a1bSJed Brown const PetscScalar QG[3] = {0.270702180178785, 0.066120127797275, -0.108402221791588};
39c4762a1bSJed Brown /* Generator constants */
40c4762a1bSJed Brown const PetscScalar H[3]    = {23.64, 6.4, 3.01};       /* Inertia constant */
41c4762a1bSJed Brown const PetscScalar Rs[3]   = {0.0, 0.0, 0.0};          /* Stator Resistance */
42c4762a1bSJed Brown const PetscScalar Xd[3]   = {0.146, 0.8958, 1.3125};  /* d-axis reactance */
43c4762a1bSJed Brown const PetscScalar Xdp[3]  = {0.0608, 0.1198, 0.1813}; /* d-axis transient reactance */
44c4762a1bSJed Brown 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 */
45c4762a1bSJed Brown const PetscScalar Xqp[3]  = {0.0969, 0.1969, 0.25};   /* q-axis transient reactance */
46c4762a1bSJed Brown const PetscScalar Td0p[3] = {8.96, 6.0, 5.89};        /* d-axis open circuit time constant */
47c4762a1bSJed Brown const PetscScalar Tq0p[3] = {0.31, 0.535, 0.6};       /* q-axis open circuit time constant */
48c4762a1bSJed Brown PetscScalar       M[3];                               /* M = 2*H/w_s */
49c4762a1bSJed Brown PetscScalar       D[3];                               /* D = 0.1*M */
50c4762a1bSJed Brown 
51c4762a1bSJed Brown PetscScalar TM[3]; /* Mechanical Torque */
52c4762a1bSJed Brown /* Exciter system constants */
53c4762a1bSJed Brown const PetscScalar KA[3] = {20.0, 20.0, 20.0};    /* Voltage regulartor gain constant */
54c4762a1bSJed Brown const PetscScalar TA[3] = {0.2, 0.2, 0.2};       /* Voltage regulator time constant */
55c4762a1bSJed Brown const PetscScalar KE[3] = {1.0, 1.0, 1.0};       /* Exciter gain constant */
56c4762a1bSJed Brown const PetscScalar TE[3] = {0.314, 0.314, 0.314}; /* Exciter time constant */
57c4762a1bSJed Brown const PetscScalar KF[3] = {0.063, 0.063, 0.063}; /* Feedback stabilizer gain constant */
58c4762a1bSJed Brown const PetscScalar TF[3] = {0.35, 0.35, 0.35};    /* Feedback stabilizer time constant */
59c4762a1bSJed Brown const PetscScalar k1[3] = {0.0039, 0.0039, 0.0039};
60c4762a1bSJed Brown const PetscScalar k2[3] = {1.555, 1.555, 1.555}; /* k1 and k2 for calculating the saturation function SE = k1*exp(k2*Efd) */
61c4762a1bSJed Brown 
62c4762a1bSJed Brown PetscScalar Vref[3];
63c4762a1bSJed Brown /* Load constants
64c4762a1bSJed Brown   We use a composite load model that describes the load and reactive powers at each time instant as follows
65c4762a1bSJed Brown   P(t) = \sum\limits_{i=0}^ld_nsegsp \ld_alphap_i*P_D0(\frac{V_m(t)}{V_m0})^\ld_betap_i
66c4762a1bSJed Brown   Q(t) = \sum\limits_{i=0}^ld_nsegsq \ld_alphaq_i*Q_D0(\frac{V_m(t)}{V_m0})^\ld_betaq_i
67c4762a1bSJed Brown   where
68c4762a1bSJed Brown     ld_nsegsp,ld_nsegsq - Number of individual load models for real and reactive power loads
69c4762a1bSJed Brown     ld_alphap,ld_alphap - Percentage contribution (weights) or loads
70c4762a1bSJed Brown     P_D0                - Real power load
71c4762a1bSJed Brown     Q_D0                - Reactive power load
72c4762a1bSJed Brown     V_m(t)              - Voltage magnitude at time t
73c4762a1bSJed Brown     V_m0                - Voltage magnitude at t = 0
74c4762a1bSJed Brown     ld_betap, ld_betaq  - exponents describing the load model for real and reactive part
75c4762a1bSJed Brown 
76c4762a1bSJed Brown     Note: All loads have the same characteristic currently.
77c4762a1bSJed Brown */
78c4762a1bSJed Brown const PetscScalar PD0[3]       = {1.25, 0.9, 1.0};
79c4762a1bSJed Brown const PetscScalar QD0[3]       = {0.5, 0.3, 0.35};
80c4762a1bSJed Brown const PetscInt    ld_nsegsp[3] = {3, 3, 3};
81c4762a1bSJed Brown const PetscScalar ld_alphap[3] = {1.0, 0.0, 0.0};
82c4762a1bSJed Brown const PetscScalar ld_betap[3]  = {2.0, 1.0, 0.0};
83c4762a1bSJed Brown const PetscInt    ld_nsegsq[3] = {3, 3, 3};
84c4762a1bSJed Brown const PetscScalar ld_alphaq[3] = {1.0, 0.0, 0.0};
85c4762a1bSJed Brown const PetscScalar ld_betaq[3]  = {2.0, 1.0, 0.0};
86c4762a1bSJed Brown 
87c4762a1bSJed Brown typedef struct {
88c4762a1bSJed Brown   DM          dmgen, dmnet;        /* DMs to manage generator and network subsystem */
89c4762a1bSJed Brown   DM          dmpgrid;             /* Composite DM to manage the entire power grid */
90c4762a1bSJed Brown   Mat         Ybus;                /* Network admittance matrix */
91c4762a1bSJed Brown   Vec         V0;                  /* Initial voltage vector (Power flow solution) */
92c4762a1bSJed Brown   PetscReal   tfaulton, tfaultoff; /* Fault on and off times */
93c4762a1bSJed Brown   PetscInt    faultbus;            /* Fault bus */
94c4762a1bSJed Brown   PetscScalar Rfault;
95c4762a1bSJed Brown   PetscReal   t0, tmax;
96c4762a1bSJed Brown   PetscInt    neqs_gen, neqs_net, neqs_pgrid;
97c4762a1bSJed Brown   Mat         Sol; /* Matrix to save solution at each time step */
98c4762a1bSJed Brown   PetscInt    stepnum;
99c4762a1bSJed Brown   PetscBool   alg_flg;
100c4762a1bSJed Brown   PetscReal   t;
101c4762a1bSJed Brown   IS          is_diff;        /* indices for differential equations */
102c4762a1bSJed Brown   IS          is_alg;         /* indices for algebraic equations */
103c4762a1bSJed Brown   PetscReal   freq_u, freq_l; /* upper and lower frequency limit */
104c4762a1bSJed Brown   PetscInt    pow;            /* power coefficient used in the cost function */
105c4762a1bSJed Brown   PetscBool   jacp_flg;
106c4762a1bSJed Brown   Mat         J, Jacp;
107c4762a1bSJed Brown   Mat         DRDU, DRDP;
108c4762a1bSJed Brown } Userctx;
109c4762a1bSJed Brown 
110c4762a1bSJed Brown /* Converts from machine frame (dq) to network (phase a real,imag) reference frame */
dq2ri(PetscScalar Fd,PetscScalar Fq,PetscScalar delta,PetscScalar * Fr,PetscScalar * Fi)111d71ae5a4SJacob Faibussowitsch PetscErrorCode dq2ri(PetscScalar Fd, PetscScalar Fq, PetscScalar delta, PetscScalar *Fr, PetscScalar *Fi)
112d71ae5a4SJacob Faibussowitsch {
113c4762a1bSJed Brown   PetscFunctionBegin;
114c4762a1bSJed Brown   *Fr = Fd * PetscSinScalar(delta) + Fq * PetscCosScalar(delta);
115c4762a1bSJed Brown   *Fi = -Fd * PetscCosScalar(delta) + Fq * PetscSinScalar(delta);
1163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
117c4762a1bSJed Brown }
118c4762a1bSJed Brown 
119c4762a1bSJed Brown /* Converts from network frame ([phase a real,imag) to machine (dq) reference frame */
ri2dq(PetscScalar Fr,PetscScalar Fi,PetscScalar delta,PetscScalar * Fd,PetscScalar * Fq)120d71ae5a4SJacob Faibussowitsch PetscErrorCode ri2dq(PetscScalar Fr, PetscScalar Fi, PetscScalar delta, PetscScalar *Fd, PetscScalar *Fq)
121d71ae5a4SJacob Faibussowitsch {
122c4762a1bSJed Brown   PetscFunctionBegin;
123c4762a1bSJed Brown   *Fd = Fr * PetscSinScalar(delta) - Fi * PetscCosScalar(delta);
124c4762a1bSJed Brown   *Fq = Fr * PetscCosScalar(delta) + Fi * PetscSinScalar(delta);
1253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
126c4762a1bSJed Brown }
127c4762a1bSJed Brown 
128c4762a1bSJed Brown /* Saves the solution at each time to a matrix */
SaveSolution(TS ts)129d71ae5a4SJacob Faibussowitsch PetscErrorCode SaveSolution(TS ts)
130d71ae5a4SJacob Faibussowitsch {
131c4762a1bSJed Brown   Userctx           *user;
132c4762a1bSJed Brown   Vec                X;
133c4762a1bSJed Brown   PetscScalar       *mat;
134c4762a1bSJed Brown   const PetscScalar *x;
135c4762a1bSJed Brown   PetscInt           idx;
136c4762a1bSJed Brown   PetscReal          t;
137c4762a1bSJed Brown 
138c4762a1bSJed Brown   PetscFunctionBegin;
1399566063dSJacob Faibussowitsch   PetscCall(TSGetApplicationContext(ts, &user));
1409566063dSJacob Faibussowitsch   PetscCall(TSGetTime(ts, &t));
1419566063dSJacob Faibussowitsch   PetscCall(TSGetSolution(ts, &X));
142c4762a1bSJed Brown   idx = user->stepnum * (user->neqs_pgrid + 1);
1439566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArray(user->Sol, &mat));
1449566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
145c4762a1bSJed Brown   mat[idx] = t;
1469566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(mat + idx + 1, x, user->neqs_pgrid));
1479566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArray(user->Sol, &mat));
1489566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
149c4762a1bSJed Brown   user->stepnum++;
1503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
151c4762a1bSJed Brown }
152c4762a1bSJed Brown 
SetInitialGuess(Vec X,Userctx * user)153d71ae5a4SJacob Faibussowitsch PetscErrorCode SetInitialGuess(Vec X, Userctx *user)
154d71ae5a4SJacob Faibussowitsch {
155c4762a1bSJed Brown   Vec          Xgen, Xnet;
156c4762a1bSJed Brown   PetscScalar *xgen, *xnet;
157c4762a1bSJed Brown   PetscInt     i, idx = 0;
158c4762a1bSJed Brown   PetscScalar  Vr, Vi, IGr, IGi, Vm, Vm2;
159c4762a1bSJed Brown   PetscScalar  Eqp, Edp, delta;
160c4762a1bSJed Brown   PetscScalar  Efd, RF, VR; /* Exciter variables */
161c4762a1bSJed Brown   PetscScalar  Id, Iq;      /* Generator dq axis currents */
162c4762a1bSJed Brown   PetscScalar  theta, Vd, Vq, SE;
163c4762a1bSJed Brown 
164c4762a1bSJed Brown   PetscFunctionBegin;
1659371c9d4SSatish Balay   M[0] = 2 * H[0] / w_s;
1669371c9d4SSatish Balay   M[1] = 2 * H[1] / w_s;
1679371c9d4SSatish Balay   M[2] = 2 * H[2] / w_s;
1689371c9d4SSatish Balay   D[0] = 0.1 * M[0];
1699371c9d4SSatish Balay   D[1] = 0.1 * M[1];
1709371c9d4SSatish Balay   D[2] = 0.1 * M[2];
171c4762a1bSJed Brown 
1729566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetLocalVectors(user->dmpgrid, &Xgen, &Xnet));
173c4762a1bSJed Brown 
174c4762a1bSJed Brown   /* Network subsystem initialization */
1759566063dSJacob Faibussowitsch   PetscCall(VecCopy(user->V0, Xnet));
176c4762a1bSJed Brown 
177c4762a1bSJed Brown   /* Generator subsystem initialization */
1789566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Xgen, &xgen));
1799566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Xnet, &xnet));
180c4762a1bSJed Brown 
181c4762a1bSJed Brown   for (i = 0; i < ngen; i++) {
182c4762a1bSJed Brown     Vr  = xnet[2 * gbus[i]];     /* Real part of generator terminal voltage */
183c4762a1bSJed Brown     Vi  = xnet[2 * gbus[i] + 1]; /* Imaginary part of the generator terminal voltage */
1849371c9d4SSatish Balay     Vm  = PetscSqrtScalar(Vr * Vr + Vi * Vi);
1859371c9d4SSatish Balay     Vm2 = Vm * Vm;
186c4762a1bSJed Brown     IGr = (Vr * PG[i] + Vi * QG[i]) / Vm2;
187c4762a1bSJed Brown     IGi = (Vi * PG[i] - Vr * QG[i]) / Vm2;
188c4762a1bSJed Brown 
189c4762a1bSJed Brown     delta = PetscAtan2Real(Vi + Xq[i] * IGr, Vr - Xq[i] * IGi); /* Machine angle */
190c4762a1bSJed Brown 
191c4762a1bSJed Brown     theta = PETSC_PI / 2.0 - delta;
192c4762a1bSJed Brown 
193c4762a1bSJed Brown     Id = IGr * PetscCosScalar(theta) - IGi * PetscSinScalar(theta); /* d-axis stator current */
194c4762a1bSJed Brown     Iq = IGr * PetscSinScalar(theta) + IGi * PetscCosScalar(theta); /* q-axis stator current */
195c4762a1bSJed Brown 
196c4762a1bSJed Brown     Vd = Vr * PetscCosScalar(theta) - Vi * PetscSinScalar(theta);
197c4762a1bSJed Brown     Vq = Vr * PetscSinScalar(theta) + Vi * PetscCosScalar(theta);
198c4762a1bSJed Brown 
199c4762a1bSJed Brown     Edp = Vd + Rs[i] * Id - Xqp[i] * Iq; /* d-axis transient EMF */
200c4762a1bSJed Brown     Eqp = Vq + Rs[i] * Iq + Xdp[i] * Id; /* q-axis transient EMF */
201c4762a1bSJed Brown 
202c4762a1bSJed Brown     TM[i] = PG[i];
203c4762a1bSJed Brown 
204c4762a1bSJed Brown     /* The generator variables are ordered as [Eqp,Edp,delta,w,Id,Iq] */
205c4762a1bSJed Brown     xgen[idx]     = Eqp;
206c4762a1bSJed Brown     xgen[idx + 1] = Edp;
207c4762a1bSJed Brown     xgen[idx + 2] = delta;
208c4762a1bSJed Brown     xgen[idx + 3] = w_s;
209c4762a1bSJed Brown 
210c4762a1bSJed Brown     idx = idx + 4;
211c4762a1bSJed Brown 
212c4762a1bSJed Brown     xgen[idx]     = Id;
213c4762a1bSJed Brown     xgen[idx + 1] = Iq;
214c4762a1bSJed Brown 
215c4762a1bSJed Brown     idx = idx + 2;
216c4762a1bSJed Brown 
217c4762a1bSJed Brown     /* Exciter */
218c4762a1bSJed Brown     Efd = Eqp + (Xd[i] - Xdp[i]) * Id;
219c4762a1bSJed Brown     SE  = k1[i] * PetscExpScalar(k2[i] * Efd);
220c4762a1bSJed Brown     VR  = KE[i] * Efd + SE;
221c4762a1bSJed Brown     RF  = KF[i] * Efd / TF[i];
222c4762a1bSJed Brown 
223c4762a1bSJed Brown     xgen[idx]     = Efd;
224c4762a1bSJed Brown     xgen[idx + 1] = RF;
225c4762a1bSJed Brown     xgen[idx + 2] = VR;
226c4762a1bSJed Brown 
227c4762a1bSJed Brown     Vref[i] = Vm + (VR / KA[i]);
228c4762a1bSJed Brown 
229c4762a1bSJed Brown     idx = idx + 3;
230c4762a1bSJed Brown   }
231c4762a1bSJed Brown 
2329566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Xgen, &xgen));
2339566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Xnet, &xnet));
234c4762a1bSJed Brown 
2359566063dSJacob Faibussowitsch   /* PetscCall(VecView(Xgen,0)); */
2369566063dSJacob Faibussowitsch   PetscCall(DMCompositeGather(user->dmpgrid, INSERT_VALUES, X, Xgen, Xnet));
2379566063dSJacob Faibussowitsch   PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid, &Xgen, &Xnet));
2383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
239c4762a1bSJed Brown }
240c4762a1bSJed Brown 
InitialGuess(Vec X,Userctx * user,const PetscScalar PGv[])241d71ae5a4SJacob Faibussowitsch PetscErrorCode InitialGuess(Vec X, Userctx *user, const PetscScalar PGv[])
242d71ae5a4SJacob Faibussowitsch {
243c4762a1bSJed Brown   Vec          Xgen, Xnet;
244c4762a1bSJed Brown   PetscScalar *xgen, *xnet;
245c4762a1bSJed Brown   PetscInt     i, idx = 0;
246c4762a1bSJed Brown   PetscScalar  Vr, Vi, IGr, IGi, Vm, Vm2;
247c4762a1bSJed Brown   PetscScalar  Eqp, Edp, delta;
248c4762a1bSJed Brown   PetscScalar  Efd, RF, VR; /* Exciter variables */
249c4762a1bSJed Brown   PetscScalar  Id, Iq;      /* Generator dq axis currents */
250c4762a1bSJed Brown   PetscScalar  theta, Vd, Vq, SE;
251c4762a1bSJed Brown 
252c4762a1bSJed Brown   PetscFunctionBegin;
2539371c9d4SSatish Balay   M[0] = 2 * H[0] / w_s;
2549371c9d4SSatish Balay   M[1] = 2 * H[1] / w_s;
2559371c9d4SSatish Balay   M[2] = 2 * H[2] / w_s;
2569371c9d4SSatish Balay   D[0] = 0.1 * M[0];
2579371c9d4SSatish Balay   D[1] = 0.1 * M[1];
2589371c9d4SSatish Balay   D[2] = 0.1 * M[2];
259c4762a1bSJed Brown 
2609566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetLocalVectors(user->dmpgrid, &Xgen, &Xnet));
261c4762a1bSJed Brown 
262c4762a1bSJed Brown   /* Network subsystem initialization */
2639566063dSJacob Faibussowitsch   PetscCall(VecCopy(user->V0, Xnet));
264c4762a1bSJed Brown 
265c4762a1bSJed Brown   /* Generator subsystem initialization */
2669566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Xgen, &xgen));
2679566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Xnet, &xnet));
268c4762a1bSJed Brown 
269c4762a1bSJed Brown   for (i = 0; i < ngen; i++) {
270c4762a1bSJed Brown     Vr  = xnet[2 * gbus[i]];     /* Real part of generator terminal voltage */
271c4762a1bSJed Brown     Vi  = xnet[2 * gbus[i] + 1]; /* Imaginary part of the generator terminal voltage */
2729371c9d4SSatish Balay     Vm  = PetscSqrtScalar(Vr * Vr + Vi * Vi);
2739371c9d4SSatish Balay     Vm2 = Vm * Vm;
274c4762a1bSJed Brown     IGr = (Vr * PGv[i] + Vi * QG[i]) / Vm2;
275c4762a1bSJed Brown     IGi = (Vi * PGv[i] - Vr * QG[i]) / Vm2;
276c4762a1bSJed Brown 
277c4762a1bSJed Brown     delta = PetscAtan2Real(Vi + Xq[i] * IGr, Vr - Xq[i] * IGi); /* Machine angle */
278c4762a1bSJed Brown 
279c4762a1bSJed Brown     theta = PETSC_PI / 2.0 - delta;
280c4762a1bSJed Brown 
281c4762a1bSJed Brown     Id = IGr * PetscCosScalar(theta) - IGi * PetscSinScalar(theta); /* d-axis stator current */
282c4762a1bSJed Brown     Iq = IGr * PetscSinScalar(theta) + IGi * PetscCosScalar(theta); /* q-axis stator current */
283c4762a1bSJed Brown 
284c4762a1bSJed Brown     Vd = Vr * PetscCosScalar(theta) - Vi * PetscSinScalar(theta);
285c4762a1bSJed Brown     Vq = Vr * PetscSinScalar(theta) + Vi * PetscCosScalar(theta);
286c4762a1bSJed Brown 
287c4762a1bSJed Brown     Edp = Vd + Rs[i] * Id - Xqp[i] * Iq; /* d-axis transient EMF */
288c4762a1bSJed Brown     Eqp = Vq + Rs[i] * Iq + Xdp[i] * Id; /* q-axis transient EMF */
289c4762a1bSJed Brown 
290c4762a1bSJed Brown     /* The generator variables are ordered as [Eqp,Edp,delta,w,Id,Iq] */
291c4762a1bSJed Brown     xgen[idx]     = Eqp;
292c4762a1bSJed Brown     xgen[idx + 1] = Edp;
293c4762a1bSJed Brown     xgen[idx + 2] = delta;
294c4762a1bSJed Brown     xgen[idx + 3] = w_s;
295c4762a1bSJed Brown 
296c4762a1bSJed Brown     idx = idx + 4;
297c4762a1bSJed Brown 
298c4762a1bSJed Brown     xgen[idx]     = Id;
299c4762a1bSJed Brown     xgen[idx + 1] = Iq;
300c4762a1bSJed Brown 
301c4762a1bSJed Brown     idx = idx + 2;
302c4762a1bSJed Brown 
303c4762a1bSJed Brown     /* Exciter */
304c4762a1bSJed Brown     Efd = Eqp + (Xd[i] - Xdp[i]) * Id;
305c4762a1bSJed Brown     SE  = k1[i] * PetscExpScalar(k2[i] * Efd);
306c4762a1bSJed Brown     VR  = KE[i] * Efd + SE;
307c4762a1bSJed Brown     RF  = KF[i] * Efd / TF[i];
308c4762a1bSJed Brown 
309c4762a1bSJed Brown     xgen[idx]     = Efd;
310c4762a1bSJed Brown     xgen[idx + 1] = RF;
311c4762a1bSJed Brown     xgen[idx + 2] = VR;
312c4762a1bSJed Brown 
313c4762a1bSJed Brown     idx = idx + 3;
314c4762a1bSJed Brown   }
315c4762a1bSJed Brown 
3169566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Xgen, &xgen));
3179566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Xnet, &xnet));
318c4762a1bSJed Brown 
3199566063dSJacob Faibussowitsch   /* PetscCall(VecView(Xgen,0)); */
3209566063dSJacob Faibussowitsch   PetscCall(DMCompositeGather(user->dmpgrid, INSERT_VALUES, X, Xgen, Xnet));
3219566063dSJacob Faibussowitsch   PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid, &Xgen, &Xnet));
3223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
323c4762a1bSJed Brown }
324c4762a1bSJed Brown 
DICDPFiniteDifference(Vec X,Vec * DICDP,Userctx * user)325d71ae5a4SJacob Faibussowitsch PetscErrorCode DICDPFiniteDifference(Vec X, Vec *DICDP, Userctx *user)
326d71ae5a4SJacob Faibussowitsch {
327c4762a1bSJed Brown   Vec         Y;
328c4762a1bSJed Brown   PetscScalar PGv[3], eps;
329c4762a1bSJed Brown   PetscInt    i, j;
330c4762a1bSJed Brown 
331362febeeSStefano Zampini   PetscFunctionBegin;
332c4762a1bSJed Brown   eps = 1.e-7;
3339566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(X, &Y));
334c4762a1bSJed Brown 
335c4762a1bSJed Brown   for (i = 0; i < ngen; i++) {
336c4762a1bSJed Brown     for (j = 0; j < 3; j++) PGv[j] = PG[j];
337c4762a1bSJed Brown     PGv[i] = PG[i] + eps;
3389566063dSJacob Faibussowitsch     PetscCall(InitialGuess(Y, user, PGv));
3399566063dSJacob Faibussowitsch     PetscCall(InitialGuess(X, user, PG));
340c4762a1bSJed Brown 
3419566063dSJacob Faibussowitsch     PetscCall(VecAXPY(Y, -1.0, X));
3429566063dSJacob Faibussowitsch     PetscCall(VecScale(Y, 1. / eps));
3439566063dSJacob Faibussowitsch     PetscCall(VecCopy(Y, DICDP[i]));
344c4762a1bSJed Brown   }
3459566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Y));
3463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
347c4762a1bSJed Brown }
348c4762a1bSJed Brown 
349c4762a1bSJed Brown /* Computes F = [-f(x,y);g(x,y)] */
ResidualFunction(SNES snes,Vec X,Vec F,Userctx * user)350d71ae5a4SJacob Faibussowitsch PetscErrorCode ResidualFunction(SNES snes, Vec X, Vec F, Userctx *user)
351d71ae5a4SJacob Faibussowitsch {
352c4762a1bSJed Brown   Vec          Xgen, Xnet, Fgen, Fnet;
353c4762a1bSJed Brown   PetscScalar *xgen, *xnet, *fgen, *fnet;
354c4762a1bSJed Brown   PetscInt     i, idx = 0;
355c4762a1bSJed Brown   PetscScalar  Vr, Vi, Vm, Vm2;
356c4762a1bSJed Brown   PetscScalar  Eqp, Edp, delta, w; /* Generator variables */
357c4762a1bSJed Brown   PetscScalar  Efd, RF, VR;        /* Exciter variables */
358c4762a1bSJed Brown   PetscScalar  Id, Iq;             /* Generator dq axis currents */
359c4762a1bSJed Brown   PetscScalar  Vd, Vq, SE;
360c4762a1bSJed Brown   PetscScalar  IGr, IGi, IDr, IDi;
361c4762a1bSJed Brown   PetscScalar  Zdq_inv[4], det;
362c4762a1bSJed Brown   PetscScalar  PD, QD, Vm0, *v0;
363c4762a1bSJed Brown   PetscInt     k;
364c4762a1bSJed Brown 
365c4762a1bSJed Brown   PetscFunctionBegin;
3669566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(F));
3679566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetLocalVectors(user->dmpgrid, &Xgen, &Xnet));
3689566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetLocalVectors(user->dmpgrid, &Fgen, &Fnet));
3699566063dSJacob Faibussowitsch   PetscCall(DMCompositeScatter(user->dmpgrid, X, Xgen, Xnet));
3709566063dSJacob Faibussowitsch   PetscCall(DMCompositeScatter(user->dmpgrid, F, Fgen, Fnet));
371c4762a1bSJed Brown 
372c4762a1bSJed Brown   /* Network current balance residual IG + Y*V + IL = 0. Only YV is added here.
373c4762a1bSJed Brown      The generator current injection, IG, and load current injection, ID are added later
374c4762a1bSJed Brown   */
375c4762a1bSJed Brown   /* Note that the values in Ybus are stored assuming the imaginary current balance
376c4762a1bSJed Brown      equation is ordered first followed by real current balance equation for each bus.
377c4762a1bSJed Brown      Thus imaginary current contribution goes in location 2*i, and
378c4762a1bSJed Brown      real current contribution in 2*i+1
379c4762a1bSJed Brown   */
3809566063dSJacob Faibussowitsch   PetscCall(MatMult(user->Ybus, Xnet, Fnet));
381c4762a1bSJed Brown 
3829566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Xgen, &xgen));
3839566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Xnet, &xnet));
3849566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Fgen, &fgen));
3859566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Fnet, &fnet));
386c4762a1bSJed Brown 
387c4762a1bSJed Brown   /* Generator subsystem */
388c4762a1bSJed Brown   for (i = 0; i < ngen; i++) {
389c4762a1bSJed Brown     Eqp   = xgen[idx];
390c4762a1bSJed Brown     Edp   = xgen[idx + 1];
391c4762a1bSJed Brown     delta = xgen[idx + 2];
392c4762a1bSJed Brown     w     = xgen[idx + 3];
393c4762a1bSJed Brown     Id    = xgen[idx + 4];
394c4762a1bSJed Brown     Iq    = xgen[idx + 5];
395c4762a1bSJed Brown     Efd   = xgen[idx + 6];
396c4762a1bSJed Brown     RF    = xgen[idx + 7];
397c4762a1bSJed Brown     VR    = xgen[idx + 8];
398c4762a1bSJed Brown 
399c4762a1bSJed Brown     /* Generator differential equations */
400c4762a1bSJed Brown     fgen[idx]     = (Eqp + (Xd[i] - Xdp[i]) * Id - Efd) / Td0p[i];
401c4762a1bSJed Brown     fgen[idx + 1] = (Edp - (Xq[i] - Xqp[i]) * Iq) / Tq0p[i];
402c4762a1bSJed Brown     fgen[idx + 2] = -w + w_s;
403c4762a1bSJed Brown     fgen[idx + 3] = (-TM[i] + Edp * Id + Eqp * Iq + (Xqp[i] - Xdp[i]) * Id * Iq + D[i] * (w - w_s)) / M[i];
404c4762a1bSJed Brown 
405c4762a1bSJed Brown     Vr = xnet[2 * gbus[i]];     /* Real part of generator terminal voltage */
406c4762a1bSJed Brown     Vi = xnet[2 * gbus[i] + 1]; /* Imaginary part of the generator terminal voltage */
407c4762a1bSJed Brown 
4089566063dSJacob Faibussowitsch     PetscCall(ri2dq(Vr, Vi, delta, &Vd, &Vq));
409c4762a1bSJed Brown     /* Algebraic equations for stator currents */
410c4762a1bSJed Brown     det = Rs[i] * Rs[i] + Xdp[i] * Xqp[i];
411c4762a1bSJed Brown 
412c4762a1bSJed Brown     Zdq_inv[0] = Rs[i] / det;
413c4762a1bSJed Brown     Zdq_inv[1] = Xqp[i] / det;
414c4762a1bSJed Brown     Zdq_inv[2] = -Xdp[i] / det;
415c4762a1bSJed Brown     Zdq_inv[3] = Rs[i] / det;
416c4762a1bSJed Brown 
417c4762a1bSJed Brown     fgen[idx + 4] = Zdq_inv[0] * (-Edp + Vd) + Zdq_inv[1] * (-Eqp + Vq) + Id;
418c4762a1bSJed Brown     fgen[idx + 5] = Zdq_inv[2] * (-Edp + Vd) + Zdq_inv[3] * (-Eqp + Vq) + Iq;
419c4762a1bSJed Brown 
420c4762a1bSJed Brown     /* Add generator current injection to network */
4219566063dSJacob Faibussowitsch     PetscCall(dq2ri(Id, Iq, delta, &IGr, &IGi));
422c4762a1bSJed Brown 
423c4762a1bSJed Brown     fnet[2 * gbus[i]] -= IGi;
424c4762a1bSJed Brown     fnet[2 * gbus[i] + 1] -= IGr;
425c4762a1bSJed Brown 
426c4762a1bSJed Brown     Vm = PetscSqrtScalar(Vd * Vd + Vq * Vq);
427c4762a1bSJed Brown 
428c4762a1bSJed Brown     SE = k1[i] * PetscExpScalar(k2[i] * Efd);
429c4762a1bSJed Brown 
430c4762a1bSJed Brown     /* Exciter differential equations */
431c4762a1bSJed Brown     fgen[idx + 6] = (KE[i] * Efd + SE - VR) / TE[i];
432c4762a1bSJed Brown     fgen[idx + 7] = (RF - KF[i] * Efd / TF[i]) / TF[i];
433c4762a1bSJed Brown     fgen[idx + 8] = (VR - KA[i] * RF + KA[i] * KF[i] * Efd / TF[i] - KA[i] * (Vref[i] - Vm)) / TA[i];
434c4762a1bSJed Brown 
435c4762a1bSJed Brown     idx = idx + 9;
436c4762a1bSJed Brown   }
437c4762a1bSJed Brown 
4389566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user->V0, &v0));
439c4762a1bSJed Brown   for (i = 0; i < nload; i++) {
440c4762a1bSJed Brown     Vr  = xnet[2 * lbus[i]];     /* Real part of load bus voltage */
441c4762a1bSJed Brown     Vi  = xnet[2 * lbus[i] + 1]; /* Imaginary part of the load bus voltage */
4429371c9d4SSatish Balay     Vm  = PetscSqrtScalar(Vr * Vr + Vi * Vi);
4439371c9d4SSatish Balay     Vm2 = Vm * Vm;
444c4762a1bSJed Brown     Vm0 = PetscSqrtScalar(v0[2 * lbus[i]] * v0[2 * lbus[i]] + v0[2 * lbus[i] + 1] * v0[2 * lbus[i] + 1]);
445c4762a1bSJed Brown     PD = QD = 0.0;
44657508eceSPierre Jolivet     for (k = 0; k < ld_nsegsp[i]; k++) PD += ld_alphap[k] * PD0[i] * PetscPowScalar(Vm / Vm0, ld_betap[k]);
44757508eceSPierre Jolivet     for (k = 0; k < ld_nsegsq[i]; k++) QD += ld_alphaq[k] * QD0[i] * PetscPowScalar(Vm / Vm0, ld_betaq[k]);
448c4762a1bSJed Brown 
449c4762a1bSJed Brown     /* Load currents */
450c4762a1bSJed Brown     IDr = (PD * Vr + QD * Vi) / Vm2;
451c4762a1bSJed Brown     IDi = (-QD * Vr + PD * Vi) / Vm2;
452c4762a1bSJed Brown 
453c4762a1bSJed Brown     fnet[2 * lbus[i]] += IDi;
454c4762a1bSJed Brown     fnet[2 * lbus[i] + 1] += IDr;
455c4762a1bSJed Brown   }
4569566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user->V0, &v0));
457c4762a1bSJed Brown 
4589566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Xgen, &xgen));
4599566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Xnet, &xnet));
4609566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Fgen, &fgen));
4619566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Fnet, &fnet));
462c4762a1bSJed Brown 
4639566063dSJacob Faibussowitsch   PetscCall(DMCompositeGather(user->dmpgrid, INSERT_VALUES, F, Fgen, Fnet));
4649566063dSJacob Faibussowitsch   PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid, &Xgen, &Xnet));
4659566063dSJacob Faibussowitsch   PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid, &Fgen, &Fnet));
4663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
467c4762a1bSJed Brown }
468c4762a1bSJed Brown 
469c4762a1bSJed Brown /* \dot{x} - f(x,y)
470c4762a1bSJed Brown      g(x,y) = 0
471c4762a1bSJed Brown  */
IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,Userctx * user)472d71ae5a4SJacob Faibussowitsch PetscErrorCode IFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, Userctx *user)
473d71ae5a4SJacob Faibussowitsch {
474c4762a1bSJed Brown   SNES               snes;
475c4762a1bSJed Brown   PetscScalar       *f;
476c4762a1bSJed Brown   const PetscScalar *xdot;
477c4762a1bSJed Brown   PetscInt           i;
478c4762a1bSJed Brown 
479c4762a1bSJed Brown   PetscFunctionBegin;
480c4762a1bSJed Brown   user->t = t;
481c4762a1bSJed Brown 
4829566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts, &snes));
4839566063dSJacob Faibussowitsch   PetscCall(ResidualFunction(snes, X, F, user));
4849566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
4859566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xdot, &xdot));
486c4762a1bSJed Brown   for (i = 0; i < ngen; i++) {
487c4762a1bSJed Brown     f[9 * i] += xdot[9 * i];
488c4762a1bSJed Brown     f[9 * i + 1] += xdot[9 * i + 1];
489c4762a1bSJed Brown     f[9 * i + 2] += xdot[9 * i + 2];
490c4762a1bSJed Brown     f[9 * i + 3] += xdot[9 * i + 3];
491c4762a1bSJed Brown     f[9 * i + 6] += xdot[9 * i + 6];
492c4762a1bSJed Brown     f[9 * i + 7] += xdot[9 * i + 7];
493c4762a1bSJed Brown     f[9 * i + 8] += xdot[9 * i + 8];
494c4762a1bSJed Brown   }
4959566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
4969566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xdot, &xdot));
4973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
498c4762a1bSJed Brown }
499c4762a1bSJed Brown 
500c4762a1bSJed Brown /* This function is used for solving the algebraic system only during fault on and
501c4762a1bSJed Brown    off times. It computes the entire F and then zeros out the part corresponding to
502c4762a1bSJed Brown    differential equations
503c4762a1bSJed Brown  F = [0;g(y)];
504c4762a1bSJed Brown */
AlgFunction(SNES snes,Vec X,Vec F,PetscCtx ctx)505*2a8381b2SBarry Smith PetscErrorCode AlgFunction(SNES snes, Vec X, Vec F, PetscCtx ctx)
506d71ae5a4SJacob Faibussowitsch {
507c4762a1bSJed Brown   Userctx     *user = (Userctx *)ctx;
508c4762a1bSJed Brown   PetscScalar *f;
509c4762a1bSJed Brown   PetscInt     i;
510c4762a1bSJed Brown 
511c4762a1bSJed Brown   PetscFunctionBegin;
5129566063dSJacob Faibussowitsch   PetscCall(ResidualFunction(snes, X, F, user));
5139566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
514c4762a1bSJed Brown   for (i = 0; i < ngen; i++) {
515c4762a1bSJed Brown     f[9 * i]     = 0;
516c4762a1bSJed Brown     f[9 * i + 1] = 0;
517c4762a1bSJed Brown     f[9 * i + 2] = 0;
518c4762a1bSJed Brown     f[9 * i + 3] = 0;
519c4762a1bSJed Brown     f[9 * i + 6] = 0;
520c4762a1bSJed Brown     f[9 * i + 7] = 0;
521c4762a1bSJed Brown     f[9 * i + 8] = 0;
522c4762a1bSJed Brown   }
5239566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
5243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
525c4762a1bSJed Brown }
526c4762a1bSJed Brown 
PreallocateJacobian(Mat J,Userctx * user)527d71ae5a4SJacob Faibussowitsch PetscErrorCode PreallocateJacobian(Mat J, Userctx *user)
528d71ae5a4SJacob Faibussowitsch {
529c4762a1bSJed Brown   PetscInt *d_nnz;
530c4762a1bSJed Brown   PetscInt  i, idx = 0, start = 0;
531c4762a1bSJed Brown   PetscInt  ncols;
532c4762a1bSJed Brown 
533c4762a1bSJed Brown   PetscFunctionBegin;
5349566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(user->neqs_pgrid, &d_nnz));
535c4762a1bSJed Brown   for (i = 0; i < user->neqs_pgrid; i++) d_nnz[i] = 0;
536c4762a1bSJed Brown   /* Generator subsystem */
537c4762a1bSJed Brown   for (i = 0; i < ngen; i++) {
538c4762a1bSJed Brown     d_nnz[idx] += 3;
539c4762a1bSJed Brown     d_nnz[idx + 1] += 2;
540c4762a1bSJed Brown     d_nnz[idx + 2] += 2;
541c4762a1bSJed Brown     d_nnz[idx + 3] += 5;
542c4762a1bSJed Brown     d_nnz[idx + 4] += 6;
543c4762a1bSJed Brown     d_nnz[idx + 5] += 6;
544c4762a1bSJed Brown 
545c4762a1bSJed Brown     d_nnz[user->neqs_gen + 2 * gbus[i]] += 3;
546c4762a1bSJed Brown     d_nnz[user->neqs_gen + 2 * gbus[i] + 1] += 3;
547c4762a1bSJed Brown 
548c4762a1bSJed Brown     d_nnz[idx + 6] += 2;
549c4762a1bSJed Brown     d_nnz[idx + 7] += 2;
550c4762a1bSJed Brown     d_nnz[idx + 8] += 5;
551c4762a1bSJed Brown 
552c4762a1bSJed Brown     idx = idx + 9;
553c4762a1bSJed Brown   }
554c4762a1bSJed Brown 
555c4762a1bSJed Brown   start = user->neqs_gen;
556c4762a1bSJed Brown   for (i = 0; i < nbus; i++) {
5579566063dSJacob Faibussowitsch     PetscCall(MatGetRow(user->Ybus, 2 * i, &ncols, NULL, NULL));
558c4762a1bSJed Brown     d_nnz[start + 2 * i] += ncols;
559c4762a1bSJed Brown     d_nnz[start + 2 * i + 1] += ncols;
5609566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(user->Ybus, 2 * i, &ncols, NULL, NULL));
561c4762a1bSJed Brown   }
562c4762a1bSJed Brown 
5639566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(J, 0, d_nnz));
5649566063dSJacob Faibussowitsch   PetscCall(PetscFree(d_nnz));
5653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
566c4762a1bSJed Brown }
567c4762a1bSJed Brown 
568c4762a1bSJed Brown /*
569c4762a1bSJed Brown    J = [-df_dx, -df_dy
570c4762a1bSJed Brown         dg_dx, dg_dy]
571c4762a1bSJed Brown */
ResidualJacobian(SNES snes,Vec X,Mat J,Mat B,PetscCtx ctx)572*2a8381b2SBarry Smith PetscErrorCode ResidualJacobian(SNES snes, Vec X, Mat J, Mat B, PetscCtx ctx)
573d71ae5a4SJacob Faibussowitsch {
574c4762a1bSJed Brown   Userctx           *user = (Userctx *)ctx;
575c4762a1bSJed Brown   Vec                Xgen, Xnet;
576c4762a1bSJed Brown   PetscScalar       *xgen, *xnet;
577c4762a1bSJed Brown   PetscInt           i, idx = 0;
578c4762a1bSJed Brown   PetscScalar        Vr, Vi, Vm, Vm2;
579c4762a1bSJed Brown   PetscScalar        Eqp, Edp, delta; /* Generator variables */
580c4762a1bSJed Brown   PetscScalar        Efd;             /* Exciter variables */
581c4762a1bSJed Brown   PetscScalar        Id, Iq;          /* Generator dq axis currents */
582c4762a1bSJed Brown   PetscScalar        Vd, Vq;
583c4762a1bSJed Brown   PetscScalar        val[10];
584c4762a1bSJed Brown   PetscInt           row[2], col[10];
585c4762a1bSJed Brown   PetscInt           net_start = user->neqs_gen;
586c4762a1bSJed Brown   PetscInt           ncols;
587c4762a1bSJed Brown   const PetscInt    *cols;
588c4762a1bSJed Brown   const PetscScalar *yvals;
589c4762a1bSJed Brown   PetscInt           k;
590c4762a1bSJed Brown   PetscScalar        Zdq_inv[4], det;
591c4762a1bSJed Brown   PetscScalar        dVd_dVr, dVd_dVi, dVq_dVr, dVq_dVi, dVd_ddelta, dVq_ddelta;
592c4762a1bSJed Brown   PetscScalar        dIGr_ddelta, dIGi_ddelta, dIGr_dId, dIGr_dIq, dIGi_dId, dIGi_dIq;
593c4762a1bSJed Brown   PetscScalar        dSE_dEfd;
594c4762a1bSJed Brown   PetscScalar        dVm_dVd, dVm_dVq, dVm_dVr, dVm_dVi;
595c4762a1bSJed Brown   PetscScalar        PD, QD, Vm0, *v0, Vm4;
596c4762a1bSJed Brown   PetscScalar        dPD_dVr, dPD_dVi, dQD_dVr, dQD_dVi;
597c4762a1bSJed Brown   PetscScalar        dIDr_dVr, dIDr_dVi, dIDi_dVr, dIDi_dVi;
598c4762a1bSJed Brown 
599c4762a1bSJed Brown   PetscFunctionBegin;
6009566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(B));
6019566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetLocalVectors(user->dmpgrid, &Xgen, &Xnet));
6029566063dSJacob Faibussowitsch   PetscCall(DMCompositeScatter(user->dmpgrid, X, Xgen, Xnet));
603c4762a1bSJed Brown 
6049566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Xgen, &xgen));
6059566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Xnet, &xnet));
606c4762a1bSJed Brown 
607c4762a1bSJed Brown   /* Generator subsystem */
608c4762a1bSJed Brown   for (i = 0; i < ngen; i++) {
609c4762a1bSJed Brown     Eqp   = xgen[idx];
610c4762a1bSJed Brown     Edp   = xgen[idx + 1];
611c4762a1bSJed Brown     delta = xgen[idx + 2];
612c4762a1bSJed Brown     Id    = xgen[idx + 4];
613c4762a1bSJed Brown     Iq    = xgen[idx + 5];
614c4762a1bSJed Brown     Efd   = xgen[idx + 6];
615c4762a1bSJed Brown 
616c4762a1bSJed Brown     /*    fgen[idx]   = (Eqp + (Xd[i] - Xdp[i])*Id - Efd)/Td0p[i]; */
617c4762a1bSJed Brown     row[0] = idx;
6189371c9d4SSatish Balay     col[0] = idx;
6199371c9d4SSatish Balay     col[1] = idx + 4;
6209371c9d4SSatish Balay     col[2] = idx + 6;
6219371c9d4SSatish Balay     val[0] = 1 / Td0p[i];
6229371c9d4SSatish Balay     val[1] = (Xd[i] - Xdp[i]) / Td0p[i];
6239371c9d4SSatish Balay     val[2] = -1 / Td0p[i];
624c4762a1bSJed Brown 
6259566063dSJacob Faibussowitsch     PetscCall(MatSetValues(J, 1, row, 3, col, val, INSERT_VALUES));
626c4762a1bSJed Brown 
627c4762a1bSJed Brown     /*    fgen[idx+1] = (Edp - (Xq[i] - Xqp[i])*Iq)/Tq0p[i]; */
628c4762a1bSJed Brown     row[0] = idx + 1;
6299371c9d4SSatish Balay     col[0] = idx + 1;
6309371c9d4SSatish Balay     col[1] = idx + 5;
6319371c9d4SSatish Balay     val[0] = 1 / Tq0p[i];
6329371c9d4SSatish Balay     val[1] = -(Xq[i] - Xqp[i]) / Tq0p[i];
6339566063dSJacob Faibussowitsch     PetscCall(MatSetValues(J, 1, row, 2, col, val, INSERT_VALUES));
634c4762a1bSJed Brown 
635c4762a1bSJed Brown     /*    fgen[idx+2] = - w + w_s; */
636c4762a1bSJed Brown     row[0] = idx + 2;
6379371c9d4SSatish Balay     col[0] = idx + 2;
6389371c9d4SSatish Balay     col[1] = idx + 3;
6399371c9d4SSatish Balay     val[0] = 0;
6409371c9d4SSatish Balay     val[1] = -1;
6419566063dSJacob Faibussowitsch     PetscCall(MatSetValues(J, 1, row, 2, col, val, INSERT_VALUES));
642c4762a1bSJed Brown 
643c4762a1bSJed Brown     /*    fgen[idx+3] = (-TM[i] + Edp*Id + Eqp*Iq + (Xqp[i] - Xdp[i])*Id*Iq + D[i]*(w - w_s))/M[i]; */
644c4762a1bSJed Brown     row[0] = idx + 3;
6459371c9d4SSatish Balay     col[0] = idx;
6469371c9d4SSatish Balay     col[1] = idx + 1;
6479371c9d4SSatish Balay     col[2] = idx + 3;
6489371c9d4SSatish Balay     col[3] = idx + 4;
6499371c9d4SSatish Balay     col[4] = idx + 5;
6509371c9d4SSatish Balay     val[0] = Iq / M[i];
6519371c9d4SSatish Balay     val[1] = Id / M[i];
6529371c9d4SSatish Balay     val[2] = D[i] / M[i];
6539371c9d4SSatish Balay     val[3] = (Edp + (Xqp[i] - Xdp[i]) * Iq) / M[i];
6549371c9d4SSatish Balay     val[4] = (Eqp + (Xqp[i] - Xdp[i]) * Id) / M[i];
6559566063dSJacob Faibussowitsch     PetscCall(MatSetValues(J, 1, row, 5, col, val, INSERT_VALUES));
656c4762a1bSJed Brown 
657c4762a1bSJed Brown     Vr = xnet[2 * gbus[i]];     /* Real part of generator terminal voltage */
658c4762a1bSJed Brown     Vi = xnet[2 * gbus[i] + 1]; /* Imaginary part of the generator terminal voltage */
6599566063dSJacob Faibussowitsch     PetscCall(ri2dq(Vr, Vi, delta, &Vd, &Vq));
660c4762a1bSJed Brown 
661c4762a1bSJed Brown     det = Rs[i] * Rs[i] + Xdp[i] * Xqp[i];
662c4762a1bSJed Brown 
663c4762a1bSJed Brown     Zdq_inv[0] = Rs[i] / det;
664c4762a1bSJed Brown     Zdq_inv[1] = Xqp[i] / det;
665c4762a1bSJed Brown     Zdq_inv[2] = -Xdp[i] / det;
666c4762a1bSJed Brown     Zdq_inv[3] = Rs[i] / det;
667c4762a1bSJed Brown 
6689371c9d4SSatish Balay     dVd_dVr    = PetscSinScalar(delta);
6699371c9d4SSatish Balay     dVd_dVi    = -PetscCosScalar(delta);
6709371c9d4SSatish Balay     dVq_dVr    = PetscCosScalar(delta);
6719371c9d4SSatish Balay     dVq_dVi    = PetscSinScalar(delta);
672c4762a1bSJed Brown     dVd_ddelta = Vr * PetscCosScalar(delta) + Vi * PetscSinScalar(delta);
673c4762a1bSJed Brown     dVq_ddelta = -Vr * PetscSinScalar(delta) + Vi * PetscCosScalar(delta);
674c4762a1bSJed Brown 
675c4762a1bSJed Brown     /*    fgen[idx+4] = Zdq_inv[0]*(-Edp + Vd) + Zdq_inv[1]*(-Eqp + Vq) + Id; */
676c4762a1bSJed Brown     row[0] = idx + 4;
6779371c9d4SSatish Balay     col[0] = idx;
6789371c9d4SSatish Balay     col[1] = idx + 1;
6799371c9d4SSatish Balay     col[2] = idx + 2;
6809371c9d4SSatish Balay     val[0] = -Zdq_inv[1];
6819371c9d4SSatish Balay     val[1] = -Zdq_inv[0];
6829371c9d4SSatish Balay     val[2] = Zdq_inv[0] * dVd_ddelta + Zdq_inv[1] * dVq_ddelta;
6839371c9d4SSatish Balay     col[3] = idx + 4;
6849371c9d4SSatish Balay     col[4] = net_start + 2 * gbus[i];
6859371c9d4SSatish Balay     col[5] = net_start + 2 * gbus[i] + 1;
6869371c9d4SSatish Balay     val[3] = 1;
6879371c9d4SSatish Balay     val[4] = Zdq_inv[0] * dVd_dVr + Zdq_inv[1] * dVq_dVr;
6889371c9d4SSatish Balay     val[5] = Zdq_inv[0] * dVd_dVi + Zdq_inv[1] * dVq_dVi;
6899566063dSJacob Faibussowitsch     PetscCall(MatSetValues(J, 1, row, 6, col, val, INSERT_VALUES));
690c4762a1bSJed Brown 
691c4762a1bSJed Brown     /*  fgen[idx+5] = Zdq_inv[2]*(-Edp + Vd) + Zdq_inv[3]*(-Eqp + Vq) + Iq; */
692c4762a1bSJed Brown     row[0] = idx + 5;
6939371c9d4SSatish Balay     col[0] = idx;
6949371c9d4SSatish Balay     col[1] = idx + 1;
6959371c9d4SSatish Balay     col[2] = idx + 2;
6969371c9d4SSatish Balay     val[0] = -Zdq_inv[3];
6979371c9d4SSatish Balay     val[1] = -Zdq_inv[2];
6989371c9d4SSatish Balay     val[2] = Zdq_inv[2] * dVd_ddelta + Zdq_inv[3] * dVq_ddelta;
6999371c9d4SSatish Balay     col[3] = idx + 5;
7009371c9d4SSatish Balay     col[4] = net_start + 2 * gbus[i];
7019371c9d4SSatish Balay     col[5] = net_start + 2 * gbus[i] + 1;
7029371c9d4SSatish Balay     val[3] = 1;
7039371c9d4SSatish Balay     val[4] = Zdq_inv[2] * dVd_dVr + Zdq_inv[3] * dVq_dVr;
7049371c9d4SSatish Balay     val[5] = Zdq_inv[2] * dVd_dVi + Zdq_inv[3] * dVq_dVi;
7059566063dSJacob Faibussowitsch     PetscCall(MatSetValues(J, 1, row, 6, col, val, INSERT_VALUES));
706c4762a1bSJed Brown 
707c4762a1bSJed Brown     dIGr_ddelta = Id * PetscCosScalar(delta) - Iq * PetscSinScalar(delta);
708c4762a1bSJed Brown     dIGi_ddelta = Id * PetscSinScalar(delta) + Iq * PetscCosScalar(delta);
7099371c9d4SSatish Balay     dIGr_dId    = PetscSinScalar(delta);
7109371c9d4SSatish Balay     dIGr_dIq    = PetscCosScalar(delta);
7119371c9d4SSatish Balay     dIGi_dId    = -PetscCosScalar(delta);
7129371c9d4SSatish Balay     dIGi_dIq    = PetscSinScalar(delta);
713c4762a1bSJed Brown 
714c4762a1bSJed Brown     /* fnet[2*gbus[i]]   -= IGi; */
715c4762a1bSJed Brown     row[0] = net_start + 2 * gbus[i];
7169371c9d4SSatish Balay     col[0] = idx + 2;
7179371c9d4SSatish Balay     col[1] = idx + 4;
7189371c9d4SSatish Balay     col[2] = idx + 5;
7199371c9d4SSatish Balay     val[0] = -dIGi_ddelta;
7209371c9d4SSatish Balay     val[1] = -dIGi_dId;
7219371c9d4SSatish Balay     val[2] = -dIGi_dIq;
7229566063dSJacob Faibussowitsch     PetscCall(MatSetValues(J, 1, row, 3, col, val, INSERT_VALUES));
723c4762a1bSJed Brown 
724c4762a1bSJed Brown     /* fnet[2*gbus[i]+1]   -= IGr; */
725c4762a1bSJed Brown     row[0] = net_start + 2 * gbus[i] + 1;
7269371c9d4SSatish Balay     col[0] = idx + 2;
7279371c9d4SSatish Balay     col[1] = idx + 4;
7289371c9d4SSatish Balay     col[2] = idx + 5;
7299371c9d4SSatish Balay     val[0] = -dIGr_ddelta;
7309371c9d4SSatish Balay     val[1] = -dIGr_dId;
7319371c9d4SSatish Balay     val[2] = -dIGr_dIq;
7329566063dSJacob Faibussowitsch     PetscCall(MatSetValues(J, 1, row, 3, col, val, INSERT_VALUES));
733c4762a1bSJed Brown 
734c4762a1bSJed Brown     Vm = PetscSqrtScalar(Vd * Vd + Vq * Vq);
735c4762a1bSJed Brown 
736c4762a1bSJed Brown     /*    fgen[idx+6] = (KE[i]*Efd + SE - VR)/TE[i]; */
737c4762a1bSJed Brown     /*    SE  = k1[i]*PetscExpScalar(k2[i]*Efd); */
738c4762a1bSJed Brown     dSE_dEfd = k1[i] * k2[i] * PetscExpScalar(k2[i] * Efd);
739c4762a1bSJed Brown 
740c4762a1bSJed Brown     row[0] = idx + 6;
7419371c9d4SSatish Balay     col[0] = idx + 6;
7429371c9d4SSatish Balay     col[1] = idx + 8;
7439371c9d4SSatish Balay     val[0] = (KE[i] + dSE_dEfd) / TE[i];
7449371c9d4SSatish Balay     val[1] = -1 / TE[i];
7459566063dSJacob Faibussowitsch     PetscCall(MatSetValues(J, 1, row, 2, col, val, INSERT_VALUES));
746c4762a1bSJed Brown 
747c4762a1bSJed Brown     /* Exciter differential equations */
748c4762a1bSJed Brown 
749c4762a1bSJed Brown     /*    fgen[idx+7] = (RF - KF[i]*Efd/TF[i])/TF[i]; */
750c4762a1bSJed Brown     row[0] = idx + 7;
7519371c9d4SSatish Balay     col[0] = idx + 6;
7529371c9d4SSatish Balay     col[1] = idx + 7;
7539371c9d4SSatish Balay     val[0] = (-KF[i] / TF[i]) / TF[i];
7549371c9d4SSatish Balay     val[1] = 1 / TF[i];
7559566063dSJacob Faibussowitsch     PetscCall(MatSetValues(J, 1, row, 2, col, val, INSERT_VALUES));
756c4762a1bSJed Brown 
757c4762a1bSJed Brown     /*    fgen[idx+8] = (VR - KA[i]*RF + KA[i]*KF[i]*Efd/TF[i] - KA[i]*(Vref[i] - Vm))/TA[i]; */
758c4762a1bSJed Brown     /* Vm = (Vd^2 + Vq^2)^0.5; */
7599371c9d4SSatish Balay     dVm_dVd = Vd / Vm;
7609371c9d4SSatish Balay     dVm_dVq = Vq / Vm;
761c4762a1bSJed Brown     dVm_dVr = dVm_dVd * dVd_dVr + dVm_dVq * dVq_dVr;
762c4762a1bSJed Brown     dVm_dVi = dVm_dVd * dVd_dVi + dVm_dVq * dVq_dVi;
763c4762a1bSJed Brown     row[0]  = idx + 8;
7649371c9d4SSatish Balay     col[0]  = idx + 6;
7659371c9d4SSatish Balay     col[1]  = idx + 7;
7669371c9d4SSatish Balay     col[2]  = idx + 8;
7679371c9d4SSatish Balay     val[0]  = (KA[i] * KF[i] / TF[i]) / TA[i];
7689371c9d4SSatish Balay     val[1]  = -KA[i] / TA[i];
7699371c9d4SSatish Balay     val[2]  = 1 / TA[i];
7709371c9d4SSatish Balay     col[3]  = net_start + 2 * gbus[i];
7719371c9d4SSatish Balay     col[4]  = net_start + 2 * gbus[i] + 1;
7729371c9d4SSatish Balay     val[3]  = KA[i] * dVm_dVr / TA[i];
7739371c9d4SSatish Balay     val[4]  = KA[i] * dVm_dVi / TA[i];
7749566063dSJacob Faibussowitsch     PetscCall(MatSetValues(J, 1, row, 5, col, val, INSERT_VALUES));
775c4762a1bSJed Brown     idx = idx + 9;
776c4762a1bSJed Brown   }
777c4762a1bSJed Brown 
778c4762a1bSJed Brown   for (i = 0; i < nbus; i++) {
7799566063dSJacob Faibussowitsch     PetscCall(MatGetRow(user->Ybus, 2 * i, &ncols, &cols, &yvals));
780c4762a1bSJed Brown     row[0] = net_start + 2 * i;
781c4762a1bSJed Brown     for (k = 0; k < ncols; k++) {
782c4762a1bSJed Brown       col[k] = net_start + cols[k];
783c4762a1bSJed Brown       val[k] = yvals[k];
784c4762a1bSJed Brown     }
7859566063dSJacob Faibussowitsch     PetscCall(MatSetValues(J, 1, row, ncols, col, val, INSERT_VALUES));
7869566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(user->Ybus, 2 * i, &ncols, &cols, &yvals));
787c4762a1bSJed Brown 
7889566063dSJacob Faibussowitsch     PetscCall(MatGetRow(user->Ybus, 2 * i + 1, &ncols, &cols, &yvals));
789c4762a1bSJed Brown     row[0] = net_start + 2 * i + 1;
790c4762a1bSJed Brown     for (k = 0; k < ncols; k++) {
791c4762a1bSJed Brown       col[k] = net_start + cols[k];
792c4762a1bSJed Brown       val[k] = yvals[k];
793c4762a1bSJed Brown     }
7949566063dSJacob Faibussowitsch     PetscCall(MatSetValues(J, 1, row, ncols, col, val, INSERT_VALUES));
7959566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(user->Ybus, 2 * i + 1, &ncols, &cols, &yvals));
796c4762a1bSJed Brown   }
797c4762a1bSJed Brown 
7989566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(J, MAT_FLUSH_ASSEMBLY));
7999566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(J, MAT_FLUSH_ASSEMBLY));
800c4762a1bSJed Brown 
8019566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user->V0, &v0));
802c4762a1bSJed Brown   for (i = 0; i < nload; i++) {
803c4762a1bSJed Brown     Vr  = xnet[2 * lbus[i]];     /* Real part of load bus voltage */
804c4762a1bSJed Brown     Vi  = xnet[2 * lbus[i] + 1]; /* Imaginary part of the load bus voltage */
8059371c9d4SSatish Balay     Vm  = PetscSqrtScalar(Vr * Vr + Vi * Vi);
8069371c9d4SSatish Balay     Vm2 = Vm * Vm;
8079371c9d4SSatish Balay     Vm4 = Vm2 * Vm2;
808c4762a1bSJed Brown     Vm0 = PetscSqrtScalar(v0[2 * lbus[i]] * v0[2 * lbus[i]] + v0[2 * lbus[i] + 1] * v0[2 * lbus[i] + 1]);
809c4762a1bSJed Brown     PD = QD = 0.0;
810c4762a1bSJed Brown     dPD_dVr = dPD_dVi = dQD_dVr = dQD_dVi = 0.0;
811c4762a1bSJed Brown     for (k = 0; k < ld_nsegsp[i]; k++) {
81257508eceSPierre Jolivet       PD += ld_alphap[k] * PD0[i] * PetscPowScalar(Vm / Vm0, ld_betap[k]);
81357508eceSPierre Jolivet       dPD_dVr += ld_alphap[k] * ld_betap[k] * PD0[i] * PetscPowScalar(1 / Vm0, ld_betap[k]) * Vr * PetscPowScalar(Vm, ld_betap[k] - 2);
81457508eceSPierre Jolivet       dPD_dVi += ld_alphap[k] * ld_betap[k] * PD0[i] * PetscPowScalar(1 / Vm0, ld_betap[k]) * Vi * PetscPowScalar(Vm, ld_betap[k] - 2);
815c4762a1bSJed Brown     }
816c4762a1bSJed Brown     for (k = 0; k < ld_nsegsq[i]; k++) {
81757508eceSPierre Jolivet       QD += ld_alphaq[k] * QD0[i] * PetscPowScalar(Vm / Vm0, ld_betaq[k]);
81857508eceSPierre Jolivet       dQD_dVr += ld_alphaq[k] * ld_betaq[k] * QD0[i] * PetscPowScalar(1 / Vm0, ld_betaq[k]) * Vr * PetscPowScalar(Vm, ld_betaq[k] - 2);
81957508eceSPierre Jolivet       dQD_dVi += ld_alphaq[k] * ld_betaq[k] * QD0[i] * PetscPowScalar(1 / Vm0, ld_betaq[k]) * Vi * PetscPowScalar(Vm, ld_betaq[k] - 2);
820c4762a1bSJed Brown     }
821c4762a1bSJed Brown 
822c4762a1bSJed Brown     /*    IDr = (PD*Vr + QD*Vi)/Vm2; */
823c4762a1bSJed Brown     /*    IDi = (-QD*Vr + PD*Vi)/Vm2; */
824c4762a1bSJed Brown 
825c4762a1bSJed Brown     dIDr_dVr = (dPD_dVr * Vr + dQD_dVr * Vi + PD) / Vm2 - ((PD * Vr + QD * Vi) * 2 * Vr) / Vm4;
826c4762a1bSJed Brown     dIDr_dVi = (dPD_dVi * Vr + dQD_dVi * Vi + QD) / Vm2 - ((PD * Vr + QD * Vi) * 2 * Vi) / Vm4;
827c4762a1bSJed Brown 
828c4762a1bSJed Brown     dIDi_dVr = (-dQD_dVr * Vr + dPD_dVr * Vi - QD) / Vm2 - ((-QD * Vr + PD * Vi) * 2 * Vr) / Vm4;
829c4762a1bSJed Brown     dIDi_dVi = (-dQD_dVi * Vr + dPD_dVi * Vi + PD) / Vm2 - ((-QD * Vr + PD * Vi) * 2 * Vi) / Vm4;
830c4762a1bSJed Brown 
831c4762a1bSJed Brown     /*    fnet[2*lbus[i]]   += IDi; */
832c4762a1bSJed Brown     row[0] = net_start + 2 * lbus[i];
8339371c9d4SSatish Balay     col[0] = net_start + 2 * lbus[i];
8349371c9d4SSatish Balay     col[1] = net_start + 2 * lbus[i] + 1;
8359371c9d4SSatish Balay     val[0] = dIDi_dVr;
8369371c9d4SSatish Balay     val[1] = dIDi_dVi;
8379566063dSJacob Faibussowitsch     PetscCall(MatSetValues(J, 1, row, 2, col, val, ADD_VALUES));
838c4762a1bSJed Brown     /*    fnet[2*lbus[i]+1] += IDr; */
839c4762a1bSJed Brown     row[0] = net_start + 2 * lbus[i] + 1;
8409371c9d4SSatish Balay     col[0] = net_start + 2 * lbus[i];
8419371c9d4SSatish Balay     col[1] = net_start + 2 * lbus[i] + 1;
8429371c9d4SSatish Balay     val[0] = dIDr_dVr;
8439371c9d4SSatish Balay     val[1] = dIDr_dVi;
8449566063dSJacob Faibussowitsch     PetscCall(MatSetValues(J, 1, row, 2, col, val, ADD_VALUES));
845c4762a1bSJed Brown   }
8469566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user->V0, &v0));
847c4762a1bSJed Brown 
8489566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Xgen, &xgen));
8499566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Xnet, &xnet));
850c4762a1bSJed Brown 
8519566063dSJacob Faibussowitsch   PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid, &Xgen, &Xnet));
852c4762a1bSJed Brown 
8539566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
8549566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
8553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
856c4762a1bSJed Brown }
857c4762a1bSJed Brown 
858c4762a1bSJed Brown /*
859c4762a1bSJed Brown    J = [I, 0
860c4762a1bSJed Brown         dg_dx, dg_dy]
861c4762a1bSJed Brown */
AlgJacobian(SNES snes,Vec X,Mat A,Mat B,PetscCtx ctx)862*2a8381b2SBarry Smith PetscErrorCode AlgJacobian(SNES snes, Vec X, Mat A, Mat B, PetscCtx ctx)
863d71ae5a4SJacob Faibussowitsch {
864c4762a1bSJed Brown   Userctx *user = (Userctx *)ctx;
865c4762a1bSJed Brown 
866c4762a1bSJed Brown   PetscFunctionBegin;
8679566063dSJacob Faibussowitsch   PetscCall(ResidualJacobian(snes, X, A, B, ctx));
8689566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A, MAT_KEEP_NONZERO_PATTERN, PETSC_TRUE));
8699566063dSJacob Faibussowitsch   PetscCall(MatZeroRowsIS(A, user->is_diff, 1.0, NULL, NULL));
8703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
871c4762a1bSJed Brown }
872c4762a1bSJed Brown 
873c4762a1bSJed Brown /*
874c4762a1bSJed Brown    J = [a*I-df_dx, -df_dy
875c4762a1bSJed Brown         dg_dx, dg_dy]
876c4762a1bSJed Brown */
877c4762a1bSJed Brown 
IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,Userctx * user)878d71ae5a4SJacob Faibussowitsch PetscErrorCode IJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat A, Mat B, Userctx *user)
879d71ae5a4SJacob Faibussowitsch {
880c4762a1bSJed Brown   SNES        snes;
881c4762a1bSJed Brown   PetscScalar atmp = (PetscScalar)a;
882c4762a1bSJed Brown   PetscInt    i, row;
883c4762a1bSJed Brown 
884c4762a1bSJed Brown   PetscFunctionBegin;
885c4762a1bSJed Brown   user->t = t;
886c4762a1bSJed Brown 
8879566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts, &snes));
8889566063dSJacob Faibussowitsch   PetscCall(ResidualJacobian(snes, X, A, B, user));
889c4762a1bSJed Brown   for (i = 0; i < ngen; i++) {
890c4762a1bSJed Brown     row = 9 * i;
8919566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES));
892c4762a1bSJed Brown     row = 9 * i + 1;
8939566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES));
894c4762a1bSJed Brown     row = 9 * i + 2;
8959566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES));
896c4762a1bSJed Brown     row = 9 * i + 3;
8979566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES));
898c4762a1bSJed Brown     row = 9 * i + 6;
8999566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES));
900c4762a1bSJed Brown     row = 9 * i + 7;
9019566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES));
902c4762a1bSJed Brown     row = 9 * i + 8;
9039566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES));
904c4762a1bSJed Brown   }
9059566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
9069566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
9073ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
908c4762a1bSJed Brown }
909c4762a1bSJed Brown 
910c4762a1bSJed Brown /* Matrix JacobianP is constant so that it only needs to be evaluated once */
RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,PetscCtx ctx0)911*2a8381b2SBarry Smith static PetscErrorCode RHSJacobianP(TS ts, PetscReal t, Vec X, Mat A, PetscCtx ctx0)
912d71ae5a4SJacob Faibussowitsch {
913c4762a1bSJed Brown   PetscScalar a;
914c4762a1bSJed Brown   PetscInt    row, col;
915c4762a1bSJed Brown   Userctx    *ctx = (Userctx *)ctx0;
916c4762a1bSJed Brown 
917c4762a1bSJed Brown   PetscFunctionBeginUser;
918c4762a1bSJed Brown   if (ctx->jacp_flg) {
9199566063dSJacob Faibussowitsch     PetscCall(MatZeroEntries(A));
920c4762a1bSJed Brown 
921c4762a1bSJed Brown     for (col = 0; col < 3; col++) {
922c4762a1bSJed Brown       a   = 1.0 / M[col];
923c4762a1bSJed Brown       row = 9 * col + 3;
9249566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A, 1, &row, 1, &col, &a, INSERT_VALUES));
925c4762a1bSJed Brown     }
926c4762a1bSJed Brown 
927c4762a1bSJed Brown     ctx->jacp_flg = PETSC_FALSE;
928c4762a1bSJed Brown 
9299566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
9309566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
931c4762a1bSJed Brown   }
9323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
933c4762a1bSJed Brown }
934c4762a1bSJed Brown 
CostIntegrand(TS ts,PetscReal t,Vec U,Vec R,Userctx * user)935d71ae5a4SJacob Faibussowitsch static PetscErrorCode CostIntegrand(TS ts, PetscReal t, Vec U, Vec R, Userctx *user)
936d71ae5a4SJacob Faibussowitsch {
937c4762a1bSJed Brown   const PetscScalar *u;
938c4762a1bSJed Brown   PetscInt           idx;
939c4762a1bSJed Brown   Vec                Xgen, Xnet;
940c4762a1bSJed Brown   PetscScalar       *r, *xgen;
941c4762a1bSJed Brown   PetscInt           i;
942c4762a1bSJed Brown 
943c4762a1bSJed Brown   PetscFunctionBegin;
9449566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetLocalVectors(user->dmpgrid, &Xgen, &Xnet));
9459566063dSJacob Faibussowitsch   PetscCall(DMCompositeScatter(user->dmpgrid, U, Xgen, Xnet));
946c4762a1bSJed Brown 
9479566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Xgen, &xgen));
948c4762a1bSJed Brown 
9499566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
9509566063dSJacob Faibussowitsch   PetscCall(VecGetArray(R, &r));
951c4762a1bSJed Brown   r[0] = 0.;
952c4762a1bSJed Brown   idx  = 0;
953c4762a1bSJed Brown   for (i = 0; i < ngen; i++) {
954c4762a1bSJed Brown     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);
955c4762a1bSJed Brown     idx += 9;
956c4762a1bSJed Brown   }
9579566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
9589566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(R, &r));
9599566063dSJacob Faibussowitsch   PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid, &Xgen, &Xnet));
9603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
961c4762a1bSJed Brown }
962c4762a1bSJed Brown 
DRDUJacobianTranspose(TS ts,PetscReal t,Vec U,Mat DRDU,Mat B,Userctx * user)963d71ae5a4SJacob Faibussowitsch static PetscErrorCode DRDUJacobianTranspose(TS ts, PetscReal t, Vec U, Mat DRDU, Mat B, Userctx *user)
964d71ae5a4SJacob Faibussowitsch {
965c4762a1bSJed Brown   Vec          Xgen, Xnet, Dgen, Dnet;
966c4762a1bSJed Brown   PetscScalar *xgen, *dgen;
967c4762a1bSJed Brown   PetscInt     i;
968c4762a1bSJed Brown   PetscInt     idx;
969c4762a1bSJed Brown   Vec          drdu_col;
970c4762a1bSJed Brown   PetscScalar *xarr;
971c4762a1bSJed Brown 
972c4762a1bSJed Brown   PetscFunctionBegin;
9739566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(U, &drdu_col));
9749566063dSJacob Faibussowitsch   PetscCall(MatDenseGetColumn(DRDU, 0, &xarr));
9759566063dSJacob Faibussowitsch   PetscCall(VecPlaceArray(drdu_col, xarr));
9769566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetLocalVectors(user->dmpgrid, &Xgen, &Xnet));
9779566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetLocalVectors(user->dmpgrid, &Dgen, &Dnet));
9789566063dSJacob Faibussowitsch   PetscCall(DMCompositeScatter(user->dmpgrid, U, Xgen, Xnet));
9799566063dSJacob Faibussowitsch   PetscCall(DMCompositeScatter(user->dmpgrid, drdu_col, Dgen, Dnet));
980c4762a1bSJed Brown 
9819566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Xgen, &xgen));
9829566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Dgen, &dgen));
983c4762a1bSJed Brown 
984c4762a1bSJed Brown   idx = 0;
985c4762a1bSJed Brown   for (i = 0; i < ngen; i++) {
986c4762a1bSJed Brown     dgen[idx + 3] = 0.;
987c4762a1bSJed Brown     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);
988c4762a1bSJed Brown     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);
989c4762a1bSJed Brown     idx += 9;
990c4762a1bSJed Brown   }
991c4762a1bSJed Brown 
9929566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Dgen, &dgen));
9939566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Xgen, &xgen));
9949566063dSJacob Faibussowitsch   PetscCall(DMCompositeGather(user->dmpgrid, INSERT_VALUES, drdu_col, Dgen, Dnet));
9959566063dSJacob Faibussowitsch   PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid, &Dgen, &Dnet));
9969566063dSJacob Faibussowitsch   PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid, &Xgen, &Xnet));
9979566063dSJacob Faibussowitsch   PetscCall(VecResetArray(drdu_col));
9989566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreColumn(DRDU, &xarr));
9999566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&drdu_col));
10003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1001c4762a1bSJed Brown }
1002c4762a1bSJed Brown 
DRDPJacobianTranspose(TS ts,PetscReal t,Vec U,Mat drdp,Userctx * user)1003d71ae5a4SJacob Faibussowitsch static PetscErrorCode DRDPJacobianTranspose(TS ts, PetscReal t, Vec U, Mat drdp, Userctx *user)
1004d71ae5a4SJacob Faibussowitsch {
1005c4762a1bSJed Brown   PetscFunctionBegin;
10063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1007c4762a1bSJed Brown }
1008c4762a1bSJed Brown 
ComputeSensiP(Vec lambda,Vec mu,Vec * DICDP,Userctx * user)1009d71ae5a4SJacob Faibussowitsch PetscErrorCode ComputeSensiP(Vec lambda, Vec mu, Vec *DICDP, Userctx *user)
1010d71ae5a4SJacob Faibussowitsch {
1011c4762a1bSJed Brown   PetscScalar *x, *y, sensip;
1012c4762a1bSJed Brown   PetscInt     i;
1013c4762a1bSJed Brown 
1014c4762a1bSJed Brown   PetscFunctionBegin;
10159566063dSJacob Faibussowitsch   PetscCall(VecGetArray(lambda, &x));
10169566063dSJacob Faibussowitsch   PetscCall(VecGetArray(mu, &y));
1017c4762a1bSJed Brown 
1018c4762a1bSJed Brown   for (i = 0; i < 3; i++) {
10199566063dSJacob Faibussowitsch     PetscCall(VecDot(lambda, DICDP[i], &sensip));
1020c4762a1bSJed Brown     sensip = sensip + y[i];
102163a3b9bcSJacob Faibussowitsch     /* PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt %" PetscInt_FMT " th parameter: %g \n",i,(double)sensip)); */
1022c4762a1bSJed Brown     y[i] = sensip;
1023c4762a1bSJed Brown   }
10249566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(mu, &y));
10253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1026c4762a1bSJed Brown }
1027c4762a1bSJed Brown 
main(int argc,char ** argv)1028d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
1029d71ae5a4SJacob Faibussowitsch {
1030c4762a1bSJed Brown   Userctx      user;
1031c4762a1bSJed Brown   Vec          p;
1032c4762a1bSJed Brown   PetscScalar *x_ptr;
1033c4762a1bSJed Brown   PetscMPIInt  size;
1034c4762a1bSJed Brown   PetscInt     i;
1035c4762a1bSJed Brown   PetscViewer  Xview, Ybusview;
1036c4762a1bSJed Brown   PetscInt    *idx2;
1037c4762a1bSJed Brown   Tao          tao;
1038c4762a1bSJed Brown   KSP          ksp;
1039c4762a1bSJed Brown   PC           pc;
1040c4762a1bSJed Brown   Vec          lowerb, upperb;
1041c4762a1bSJed Brown 
1042327415f7SBarry Smith   PetscFunctionBeginUser;
10439566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, "petscoptions", help));
10449566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
10453c633725SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs");
1046c4762a1bSJed Brown 
1047c4762a1bSJed Brown   user.jacp_flg   = PETSC_TRUE;
1048c4762a1bSJed Brown   user.neqs_gen   = 9 * ngen; /* # eqs. for generator subsystem */
1049c4762a1bSJed Brown   user.neqs_net   = 2 * nbus; /* # eqs. for network subsystem   */
1050c4762a1bSJed Brown   user.neqs_pgrid = user.neqs_gen + user.neqs_net;
1051c4762a1bSJed Brown 
1052c4762a1bSJed Brown   /* Create indices for differential and algebraic equations */
10539566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(7 * ngen, &idx2));
1054c4762a1bSJed Brown   for (i = 0; i < ngen; i++) {
10559371c9d4SSatish Balay     idx2[7 * i]     = 9 * i;
10569371c9d4SSatish Balay     idx2[7 * i + 1] = 9 * i + 1;
10579371c9d4SSatish Balay     idx2[7 * i + 2] = 9 * i + 2;
10589371c9d4SSatish Balay     idx2[7 * i + 3] = 9 * i + 3;
10599371c9d4SSatish Balay     idx2[7 * i + 4] = 9 * i + 6;
10609371c9d4SSatish Balay     idx2[7 * i + 5] = 9 * i + 7;
10619371c9d4SSatish Balay     idx2[7 * i + 6] = 9 * i + 8;
1062c4762a1bSJed Brown   }
10639566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, 7 * ngen, idx2, PETSC_COPY_VALUES, &user.is_diff));
10649566063dSJacob Faibussowitsch   PetscCall(ISComplement(user.is_diff, 0, user.neqs_pgrid, &user.is_alg));
10659566063dSJacob Faibussowitsch   PetscCall(PetscFree(idx2));
1066c4762a1bSJed Brown 
1067c4762a1bSJed Brown   /* Set run time options */
1068d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Transient stability fault options", "");
1069c4762a1bSJed Brown   {
1070c4762a1bSJed Brown     user.tfaulton  = 1.0;
1071c4762a1bSJed Brown     user.tfaultoff = 1.2;
1072c4762a1bSJed Brown     user.Rfault    = 0.0001;
1073c4762a1bSJed Brown     user.faultbus  = 8;
10749566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-tfaulton", "", "", user.tfaulton, &user.tfaulton, NULL));
10759566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-tfaultoff", "", "", user.tfaultoff, &user.tfaultoff, NULL));
10769566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-faultbus", "", "", user.faultbus, &user.faultbus, NULL));
1077c4762a1bSJed Brown     user.t0   = 0.0;
1078c4762a1bSJed Brown     user.tmax = 1.3;
10799566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-t0", "", "", user.t0, &user.t0, NULL));
10809566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-tmax", "", "", user.tmax, &user.tmax, NULL));
1081c4762a1bSJed Brown     user.freq_u = 61.0;
1082c4762a1bSJed Brown     user.freq_l = 59.0;
1083c4762a1bSJed Brown     user.pow    = 2;
10849566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-frequ", "", "", user.freq_u, &user.freq_u, NULL));
10859566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-freql", "", "", user.freq_l, &user.freq_l, NULL));
10869566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-pow", "", "", user.pow, &user.pow, NULL));
1087c4762a1bSJed Brown   }
1088d0609cedSBarry Smith   PetscOptionsEnd();
1089c4762a1bSJed Brown 
1090c4762a1bSJed Brown   /* Create DMs for generator and network subsystems */
10919566063dSJacob Faibussowitsch   PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, user.neqs_gen, 1, 1, NULL, &user.dmgen));
10929566063dSJacob Faibussowitsch   PetscCall(DMSetOptionsPrefix(user.dmgen, "dmgen_"));
10939566063dSJacob Faibussowitsch   PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, user.neqs_net, 1, 1, NULL, &user.dmnet));
10949566063dSJacob Faibussowitsch   PetscCall(DMSetOptionsPrefix(user.dmnet, "dmnet_"));
10959566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(user.dmnet));
10969566063dSJacob Faibussowitsch   PetscCall(DMSetUp(user.dmnet));
1097c4762a1bSJed Brown   /* Create a composite DM packer and add the two DMs */
10989566063dSJacob Faibussowitsch   PetscCall(DMCompositeCreate(PETSC_COMM_WORLD, &user.dmpgrid));
10999566063dSJacob Faibussowitsch   PetscCall(DMSetOptionsPrefix(user.dmpgrid, "pgrid_"));
11009566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(user.dmgen));
11019566063dSJacob Faibussowitsch   PetscCall(DMSetUp(user.dmgen));
11029566063dSJacob Faibussowitsch   PetscCall(DMCompositeAddDM(user.dmpgrid, user.dmgen));
11039566063dSJacob Faibussowitsch   PetscCall(DMCompositeAddDM(user.dmpgrid, user.dmnet));
1104c4762a1bSJed Brown 
1105c4762a1bSJed Brown   /* Read initial voltage vector and Ybus */
11069566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "X.bin", FILE_MODE_READ, &Xview));
11079566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "Ybus.bin", FILE_MODE_READ, &Ybusview));
1108c4762a1bSJed Brown 
11099566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &user.V0));
11109566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(user.V0, PETSC_DECIDE, user.neqs_net));
11119566063dSJacob Faibussowitsch   PetscCall(VecLoad(user.V0, Xview));
1112c4762a1bSJed Brown 
11139566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &user.Ybus));
11149566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(user.Ybus, PETSC_DECIDE, PETSC_DECIDE, user.neqs_net, user.neqs_net));
11159566063dSJacob Faibussowitsch   PetscCall(MatSetType(user.Ybus, MATBAIJ));
11169566063dSJacob Faibussowitsch   /*  PetscCall(MatSetBlockSize(ctx->Ybus,2)); */
11179566063dSJacob Faibussowitsch   PetscCall(MatLoad(user.Ybus, Ybusview));
1118c4762a1bSJed Brown 
11199566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&Xview));
11209566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&Ybusview));
1121c4762a1bSJed Brown 
1122c4762a1bSJed Brown   /* Allocate space for Jacobians */
11239566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &user.J));
11249566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(user.J, PETSC_DECIDE, PETSC_DECIDE, user.neqs_pgrid, user.neqs_pgrid));
11259566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(user.J));
11269566063dSJacob Faibussowitsch   PetscCall(PreallocateJacobian(user.J, &user));
1127c4762a1bSJed Brown 
11289566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &user.Jacp));
11299566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(user.Jacp, PETSC_DECIDE, PETSC_DECIDE, user.neqs_pgrid, 3));
11309566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(user.Jacp));
11319566063dSJacob Faibussowitsch   PetscCall(MatSetUp(user.Jacp));
11329566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(user.Jacp)); /* initialize to zeros */
1133c4762a1bSJed Brown 
11349566063dSJacob Faibussowitsch   PetscCall(MatCreateDense(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, 3, 1, NULL, &user.DRDP));
11359566063dSJacob Faibussowitsch   PetscCall(MatSetUp(user.DRDP));
11369566063dSJacob Faibussowitsch   PetscCall(MatCreateDense(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, user.neqs_pgrid, 1, NULL, &user.DRDU));
11379566063dSJacob Faibussowitsch   PetscCall(MatSetUp(user.DRDU));
1138c4762a1bSJed Brown 
1139c4762a1bSJed Brown   /* Create TAO solver and set desired solution method */
11409566063dSJacob Faibussowitsch   PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao));
11419566063dSJacob Faibussowitsch   PetscCall(TaoSetType(tao, TAOBLMVM));
1142c4762a1bSJed Brown   /*
1143c4762a1bSJed Brown      Optimization starts
1144c4762a1bSJed Brown   */
1145c4762a1bSJed Brown   /* Set initial solution guess */
11469566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_WORLD, 3, &p));
11479566063dSJacob Faibussowitsch   PetscCall(VecGetArray(p, &x_ptr));
11489371c9d4SSatish Balay   x_ptr[0] = PG[0];
11499371c9d4SSatish Balay   x_ptr[1] = PG[1];
11509371c9d4SSatish Balay   x_ptr[2] = PG[2];
11519566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(p, &x_ptr));
1152c4762a1bSJed Brown 
11539566063dSJacob Faibussowitsch   PetscCall(TaoSetSolution(tao, p));
1154c4762a1bSJed Brown   /* Set routine for function and gradient evaluation */
11559566063dSJacob Faibussowitsch   PetscCall(TaoSetObjectiveAndGradient(tao, NULL, FormFunctionGradient, &user));
1156c4762a1bSJed Brown 
1157c4762a1bSJed Brown   /* Set bounds for the optimization */
11589566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(p, &lowerb));
11599566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(p, &upperb));
11609566063dSJacob Faibussowitsch   PetscCall(VecGetArray(lowerb, &x_ptr));
11619371c9d4SSatish Balay   x_ptr[0] = 0.5;
11629371c9d4SSatish Balay   x_ptr[1] = 0.5;
11639371c9d4SSatish Balay   x_ptr[2] = 0.5;
11649566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(lowerb, &x_ptr));
11659566063dSJacob Faibussowitsch   PetscCall(VecGetArray(upperb, &x_ptr));
11669371c9d4SSatish Balay   x_ptr[0] = 2.0;
11679371c9d4SSatish Balay   x_ptr[1] = 2.0;
11689371c9d4SSatish Balay   x_ptr[2] = 2.0;
11699566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(upperb, &x_ptr));
11709566063dSJacob Faibussowitsch   PetscCall(TaoSetVariableBounds(tao, lowerb, upperb));
1171c4762a1bSJed Brown 
1172c4762a1bSJed Brown   /* Check for any TAO command line options */
11739566063dSJacob Faibussowitsch   PetscCall(TaoSetFromOptions(tao));
11749566063dSJacob Faibussowitsch   PetscCall(TaoGetKSP(tao, &ksp));
1175c4762a1bSJed Brown   if (ksp) {
11769566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(ksp, &pc));
11779566063dSJacob Faibussowitsch     PetscCall(PCSetType(pc, PCNONE));
1178c4762a1bSJed Brown   }
1179c4762a1bSJed Brown 
1180c4762a1bSJed Brown   /* SOLVE THE APPLICATION */
11819566063dSJacob Faibussowitsch   PetscCall(TaoSolve(tao));
1182c4762a1bSJed Brown 
11839566063dSJacob Faibussowitsch   PetscCall(VecView(p, PETSC_VIEWER_STDOUT_WORLD));
1184c4762a1bSJed Brown   /* Free TAO data structures */
11859566063dSJacob Faibussowitsch   PetscCall(TaoDestroy(&tao));
1186c4762a1bSJed Brown 
11879566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&user.dmgen));
11889566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&user.dmnet));
11899566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&user.dmpgrid));
11909566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&user.is_diff));
11919566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&user.is_alg));
1192c4762a1bSJed Brown 
11939566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user.J));
11949566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user.Jacp));
11959566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user.Ybus));
11969566063dSJacob Faibussowitsch   /* PetscCall(MatDestroy(&user.Sol)); */
11979566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.V0));
11989566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&p));
11999566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&lowerb));
12009566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&upperb));
12019566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user.DRDU));
12029566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user.DRDP));
12039566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
1204b122ec5aSJacob Faibussowitsch   return 0;
1205c4762a1bSJed Brown }
1206c4762a1bSJed Brown 
1207c4762a1bSJed Brown /* ------------------------------------------------------------------ */
1208c4762a1bSJed Brown /*
1209c4762a1bSJed Brown    FormFunction - Evaluates the function and corresponding gradient.
1210c4762a1bSJed Brown 
1211c4762a1bSJed Brown    Input Parameters:
1212c4762a1bSJed Brown    tao - the Tao context
1213c4762a1bSJed Brown    X   - the input vector
1214a82e8c82SStefano Zampini    ptr - optional user-defined context, as set by TaoSetObjectiveAndGradient()
1215c4762a1bSJed Brown 
1216c4762a1bSJed Brown    Output Parameters:
1217c4762a1bSJed Brown    f   - the newly evaluated function
1218c4762a1bSJed Brown    G   - the newly evaluated gradient
1219c4762a1bSJed Brown */
FormFunctionGradient(Tao tao,Vec P,PetscReal * f,Vec G,PetscCtx ctx0)1220*2a8381b2SBarry Smith PetscErrorCode FormFunctionGradient(Tao tao, Vec P, PetscReal *f, Vec G, PetscCtx ctx0)
1221d71ae5a4SJacob Faibussowitsch {
1222c4762a1bSJed Brown   TS       ts, quadts;
1223c4762a1bSJed Brown   SNES     snes_alg;
1224c4762a1bSJed Brown   Userctx *ctx = (Userctx *)ctx0;
1225c4762a1bSJed Brown   Vec      X;
1226c4762a1bSJed Brown   PetscInt i;
1227c4762a1bSJed Brown   /* sensitivity context */
1228c4762a1bSJed Brown   PetscScalar *x_ptr;
1229c4762a1bSJed Brown   Vec          lambda[1], q;
1230c4762a1bSJed Brown   Vec          mu[1];
1231c4762a1bSJed Brown   PetscInt     steps1, steps2, steps3;
1232c4762a1bSJed Brown   Vec          DICDP[3];
1233c4762a1bSJed Brown   Vec          F_alg;
1234c4762a1bSJed Brown   PetscInt     row_loc, col_loc;
1235c4762a1bSJed Brown   PetscScalar  val;
1236c4762a1bSJed Brown   Vec          Xdot;
1237c4762a1bSJed Brown 
1238c4762a1bSJed Brown   PetscFunctionBegin;
12399566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(P, (const PetscScalar **)&x_ptr));
1240c4762a1bSJed Brown   PG[0] = x_ptr[0];
1241c4762a1bSJed Brown   PG[1] = x_ptr[1];
1242c4762a1bSJed Brown   PG[2] = x_ptr[2];
12439566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(P, (const PetscScalar **)&x_ptr));
1244c4762a1bSJed Brown 
1245c4762a1bSJed Brown   ctx->stepnum = 0;
1246c4762a1bSJed Brown 
12479566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(ctx->dmpgrid, &X));
1248c4762a1bSJed Brown 
1249c4762a1bSJed Brown   /* Create matrix to save solutions at each time step */
12509566063dSJacob Faibussowitsch   /* PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,ctx->neqs_pgrid+1,1002,NULL,&ctx->Sol)); */
1251c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1252c4762a1bSJed Brown      Create timestepping solver context
1253c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
12549566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
12559566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
12569566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSCN));
12578434afd1SBarry Smith   PetscCall(TSSetIFunction(ts, NULL, (TSIFunctionFn *)IFunction, ctx));
12588434afd1SBarry Smith   PetscCall(TSSetIJacobian(ts, ctx->J, ctx->J, (TSIJacobianFn *)IJacobian, ctx));
12599566063dSJacob Faibussowitsch   PetscCall(TSSetApplicationContext(ts, ctx));
1260c4762a1bSJed Brown   /*   Set RHS JacobianP */
12619566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobianP(ts, ctx->Jacp, RHSJacobianP, ctx));
1262c4762a1bSJed Brown 
12639566063dSJacob Faibussowitsch   PetscCall(TSCreateQuadratureTS(ts, PETSC_FALSE, &quadts));
12648434afd1SBarry Smith   PetscCall(TSSetRHSFunction(quadts, NULL, (TSRHSFunctionFn *)CostIntegrand, ctx));
12658434afd1SBarry Smith   PetscCall(TSSetRHSJacobian(quadts, ctx->DRDU, ctx->DRDU, (TSRHSJacobianFn *)DRDUJacobianTranspose, ctx));
12668434afd1SBarry Smith   PetscCall(TSSetRHSJacobianP(quadts, ctx->DRDP, (TSRHSJacobianPFn *)DRDPJacobianTranspose, ctx));
1267c4762a1bSJed Brown 
1268c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1269c4762a1bSJed Brown      Set initial conditions
1270c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
12719566063dSJacob Faibussowitsch   PetscCall(SetInitialGuess(X, ctx));
1272c4762a1bSJed Brown 
1273c4762a1bSJed Brown   /* Approximate DICDP with finite difference, we want to zero out network variables */
127448a46eb9SPierre Jolivet   for (i = 0; i < 3; i++) PetscCall(VecDuplicate(X, &DICDP[i]));
12759566063dSJacob Faibussowitsch   PetscCall(DICDPFiniteDifference(X, DICDP, ctx));
1276c4762a1bSJed Brown 
12779566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(X, &F_alg));
12789566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes_alg));
12799566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes_alg, F_alg, AlgFunction, ctx));
12809566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(ctx->J));
12819566063dSJacob Faibussowitsch   PetscCall(SNESSetJacobian(snes_alg, ctx->J, ctx->J, AlgJacobian, ctx));
12829566063dSJacob Faibussowitsch   PetscCall(SNESSetOptionsPrefix(snes_alg, "alg_"));
12839566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes_alg));
1284c4762a1bSJed Brown   ctx->alg_flg = PETSC_TRUE;
1285c4762a1bSJed Brown   /* Solve the algebraic equations */
12869566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes_alg, NULL, X));
1287c4762a1bSJed Brown 
1288c4762a1bSJed Brown   /* Just to set up the Jacobian structure */
12899566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(X, &Xdot));
12909566063dSJacob Faibussowitsch   PetscCall(IJacobian(ts, 0.0, X, Xdot, 0.0, ctx->J, ctx->J, ctx));
12919566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Xdot));
1292c4762a1bSJed Brown 
1293c4762a1bSJed Brown   ctx->stepnum++;
1294c4762a1bSJed Brown 
1295c4762a1bSJed Brown   /*
1296c4762a1bSJed Brown     Save trajectory of solution so that TSAdjointSolve() may be used
1297c4762a1bSJed Brown   */
12989566063dSJacob Faibussowitsch   PetscCall(TSSetSaveTrajectory(ts));
1299c4762a1bSJed Brown 
13009566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, 0.01));
13019566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
13029566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
13039566063dSJacob Faibussowitsch   /* PetscCall(TSSetPostStep(ts,SaveSolution)); */
1304c4762a1bSJed Brown 
1305c4762a1bSJed Brown   /* Prefault period */
1306c4762a1bSJed Brown   ctx->alg_flg = PETSC_FALSE;
13079566063dSJacob Faibussowitsch   PetscCall(TSSetTime(ts, 0.0));
13089566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, ctx->tfaulton));
13099566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, X));
13109566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts, &steps1));
1311c4762a1bSJed Brown 
1312c4762a1bSJed Brown   /* Create the nonlinear solver for solving the algebraic system */
1313c4762a1bSJed Brown   /* Note that although the algebraic system needs to be solved only for
1314c4762a1bSJed Brown      Idq and V, we reuse the entire system including xgen. The xgen
1315c4762a1bSJed Brown      variables are held constant by setting their residuals to 0 and
1316c4762a1bSJed Brown      putting a 1 on the Jacobian diagonal for xgen rows
1317c4762a1bSJed Brown   */
13189566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(ctx->J));
1319c4762a1bSJed Brown 
1320c4762a1bSJed Brown   /* Apply disturbance - resistive fault at ctx->faultbus */
1321c4762a1bSJed Brown   /* This is done by adding shunt conductance to the diagonal location
1322c4762a1bSJed Brown      in the Ybus matrix */
13239371c9d4SSatish Balay   row_loc = 2 * ctx->faultbus;
13249371c9d4SSatish Balay   col_loc = 2 * ctx->faultbus + 1; /* Location for G */
1325c4762a1bSJed Brown   val     = 1 / ctx->Rfault;
13269566063dSJacob Faibussowitsch   PetscCall(MatSetValues(ctx->Ybus, 1, &row_loc, 1, &col_loc, &val, ADD_VALUES));
13279371c9d4SSatish Balay   row_loc = 2 * ctx->faultbus + 1;
13289371c9d4SSatish Balay   col_loc = 2 * ctx->faultbus; /* Location for G */
1329c4762a1bSJed Brown   val     = 1 / ctx->Rfault;
13309566063dSJacob Faibussowitsch   PetscCall(MatSetValues(ctx->Ybus, 1, &row_loc, 1, &col_loc, &val, ADD_VALUES));
1331c4762a1bSJed Brown 
13329566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(ctx->Ybus, MAT_FINAL_ASSEMBLY));
13339566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(ctx->Ybus, MAT_FINAL_ASSEMBLY));
1334c4762a1bSJed Brown 
1335c4762a1bSJed Brown   ctx->alg_flg = PETSC_TRUE;
1336c4762a1bSJed Brown   /* Solve the algebraic equations */
13379566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes_alg, NULL, X));
1338c4762a1bSJed Brown 
1339c4762a1bSJed Brown   ctx->stepnum++;
1340c4762a1bSJed Brown 
1341c4762a1bSJed Brown   /* Disturbance period */
1342c4762a1bSJed Brown   ctx->alg_flg = PETSC_FALSE;
13439566063dSJacob Faibussowitsch   PetscCall(TSSetTime(ts, ctx->tfaulton));
13449566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, ctx->tfaultoff));
13459566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, X));
13469566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts, &steps2));
1347c4762a1bSJed Brown   steps2 -= steps1;
1348c4762a1bSJed Brown 
1349c4762a1bSJed Brown   /* Remove the fault */
13509371c9d4SSatish Balay   row_loc = 2 * ctx->faultbus;
13519371c9d4SSatish Balay   col_loc = 2 * ctx->faultbus + 1;
1352c4762a1bSJed Brown   val     = -1 / ctx->Rfault;
13539566063dSJacob Faibussowitsch   PetscCall(MatSetValues(ctx->Ybus, 1, &row_loc, 1, &col_loc, &val, ADD_VALUES));
13549371c9d4SSatish Balay   row_loc = 2 * ctx->faultbus + 1;
13559371c9d4SSatish Balay   col_loc = 2 * ctx->faultbus;
1356c4762a1bSJed Brown   val     = -1 / ctx->Rfault;
13579566063dSJacob Faibussowitsch   PetscCall(MatSetValues(ctx->Ybus, 1, &row_loc, 1, &col_loc, &val, ADD_VALUES));
1358c4762a1bSJed Brown 
13599566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(ctx->Ybus, MAT_FINAL_ASSEMBLY));
13609566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(ctx->Ybus, MAT_FINAL_ASSEMBLY));
1361c4762a1bSJed Brown 
13629566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(ctx->J));
1363c4762a1bSJed Brown 
1364c4762a1bSJed Brown   ctx->alg_flg = PETSC_TRUE;
1365c4762a1bSJed Brown 
1366c4762a1bSJed Brown   /* Solve the algebraic equations */
13679566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes_alg, NULL, X));
1368c4762a1bSJed Brown 
1369c4762a1bSJed Brown   ctx->stepnum++;
1370c4762a1bSJed Brown 
1371c4762a1bSJed Brown   /* Post-disturbance period */
1372c4762a1bSJed Brown   ctx->alg_flg = PETSC_TRUE;
13739566063dSJacob Faibussowitsch   PetscCall(TSSetTime(ts, ctx->tfaultoff));
13749566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, ctx->tmax));
13759566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, X));
13769566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts, &steps3));
1377c4762a1bSJed Brown   steps3 -= steps2;
1378c4762a1bSJed Brown   steps3 -= steps1;
1379c4762a1bSJed Brown 
1380c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1381c4762a1bSJed Brown      Adjoint model starts here
1382c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
13839566063dSJacob Faibussowitsch   PetscCall(TSSetPostStep(ts, NULL));
13849566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(ctx->J, &lambda[0], NULL));
1385c4762a1bSJed Brown   /*   Set initial conditions for the adjoint integration */
13869566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(lambda[0]));
1387c4762a1bSJed Brown 
13889566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(ctx->Jacp, &mu[0], NULL));
13899566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(mu[0]));
13909566063dSJacob Faibussowitsch   PetscCall(TSSetCostGradients(ts, 1, lambda, mu));
1391c4762a1bSJed Brown 
13929566063dSJacob Faibussowitsch   PetscCall(TSAdjointSetSteps(ts, steps3));
13939566063dSJacob Faibussowitsch   PetscCall(TSAdjointSolve(ts));
1394c4762a1bSJed Brown 
13959566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(ctx->J));
1396c4762a1bSJed Brown   /* Applying disturbance - resistive fault at ctx->faultbus */
1397c4762a1bSJed Brown   /* This is done by deducting shunt conductance to the diagonal location
1398c4762a1bSJed Brown      in the Ybus matrix */
13999371c9d4SSatish Balay   row_loc = 2 * ctx->faultbus;
14009371c9d4SSatish Balay   col_loc = 2 * ctx->faultbus + 1; /* Location for G */
1401c4762a1bSJed Brown   val     = 1. / ctx->Rfault;
14029566063dSJacob Faibussowitsch   PetscCall(MatSetValues(ctx->Ybus, 1, &row_loc, 1, &col_loc, &val, ADD_VALUES));
14039371c9d4SSatish Balay   row_loc = 2 * ctx->faultbus + 1;
14049371c9d4SSatish Balay   col_loc = 2 * ctx->faultbus; /* Location for G */
1405c4762a1bSJed Brown   val     = 1. / ctx->Rfault;
14069566063dSJacob Faibussowitsch   PetscCall(MatSetValues(ctx->Ybus, 1, &row_loc, 1, &col_loc, &val, ADD_VALUES));
1407c4762a1bSJed Brown 
14089566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(ctx->Ybus, MAT_FINAL_ASSEMBLY));
14099566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(ctx->Ybus, MAT_FINAL_ASSEMBLY));
1410c4762a1bSJed Brown 
1411c4762a1bSJed Brown   /*   Set number of steps for the adjoint integration */
14129566063dSJacob Faibussowitsch   PetscCall(TSAdjointSetSteps(ts, steps2));
14139566063dSJacob Faibussowitsch   PetscCall(TSAdjointSolve(ts));
1414c4762a1bSJed Brown 
14159566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(ctx->J));
1416c4762a1bSJed Brown   /* remove the fault */
14179371c9d4SSatish Balay   row_loc = 2 * ctx->faultbus;
14189371c9d4SSatish Balay   col_loc = 2 * ctx->faultbus + 1; /* Location for G */
1419c4762a1bSJed Brown   val     = -1. / ctx->Rfault;
14209566063dSJacob Faibussowitsch   PetscCall(MatSetValues(ctx->Ybus, 1, &row_loc, 1, &col_loc, &val, ADD_VALUES));
14219371c9d4SSatish Balay   row_loc = 2 * ctx->faultbus + 1;
14229371c9d4SSatish Balay   col_loc = 2 * ctx->faultbus; /* Location for G */
1423c4762a1bSJed Brown   val     = -1. / ctx->Rfault;
14249566063dSJacob Faibussowitsch   PetscCall(MatSetValues(ctx->Ybus, 1, &row_loc, 1, &col_loc, &val, ADD_VALUES));
1425c4762a1bSJed Brown 
14269566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(ctx->Ybus, MAT_FINAL_ASSEMBLY));
14279566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(ctx->Ybus, MAT_FINAL_ASSEMBLY));
1428c4762a1bSJed Brown 
1429c4762a1bSJed Brown   /*   Set number of steps for the adjoint integration */
14309566063dSJacob Faibussowitsch   PetscCall(TSAdjointSetSteps(ts, steps1));
14319566063dSJacob Faibussowitsch   PetscCall(TSAdjointSolve(ts));
1432c4762a1bSJed Brown 
14339566063dSJacob Faibussowitsch   PetscCall(ComputeSensiP(lambda[0], mu[0], DICDP, ctx));
14349566063dSJacob Faibussowitsch   PetscCall(VecCopy(mu[0], G));
1435c4762a1bSJed Brown 
14369566063dSJacob Faibussowitsch   PetscCall(TSGetQuadratureTS(ts, NULL, &quadts));
14379566063dSJacob Faibussowitsch   PetscCall(TSGetSolution(quadts, &q));
14389566063dSJacob Faibussowitsch   PetscCall(VecGetArray(q, &x_ptr));
1439c4762a1bSJed Brown   *f       = x_ptr[0];
1440c4762a1bSJed Brown   x_ptr[0] = 0;
14419566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(q, &x_ptr));
1442c4762a1bSJed Brown 
14439566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&lambda[0]));
14449566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mu[0]));
1445c4762a1bSJed Brown 
14469566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes_alg));
14479566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&F_alg));
14489566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X));
14499566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
145048a46eb9SPierre Jolivet   for (i = 0; i < 3; i++) PetscCall(VecDestroy(&DICDP[i]));
14513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1452c4762a1bSJed Brown }
1453c4762a1bSJed Brown 
1454c4762a1bSJed Brown /*TEST
1455c4762a1bSJed Brown 
1456c4762a1bSJed Brown    build:
1457dfd57a17SPierre Jolivet       requires: double !complex !defined(PETSC_USE_64BIT_INDICES)
1458c4762a1bSJed Brown 
1459c4762a1bSJed Brown    test:
1460c4762a1bSJed Brown       args: -viewer_binary_skip_info -tao_monitor -tao_gttol .2
1461c4762a1bSJed Brown       localrunfiles: petscoptions X.bin Ybus.bin
1462c4762a1bSJed Brown 
1463c4762a1bSJed Brown TEST*/
1464