1c4762a1bSJed Brown static char help[] = "Sensitivity analysis applied in 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 The equations for the stability analysis are described by the DAE. See ex9bus.c for details.
10c4762a1bSJed Brown The system has 'jumps' due to faults, thus the time interval is split into multiple sections, and TSSolve is called for each of them. But TSAdjointSolve only needs to be called once since the whole trajectory has been saved in the forward run.
11c4762a1bSJed Brown The code computes the sensitivity of a final state w.r.t. initial conditions.
12c4762a1bSJed Brown */
13c4762a1bSJed Brown
14c4762a1bSJed Brown #include <petscts.h>
15c4762a1bSJed Brown #include <petscdm.h>
16c4762a1bSJed Brown #include <petscdmda.h>
17c4762a1bSJed Brown #include <petscdmcomposite.h>
18c4762a1bSJed Brown
19c4762a1bSJed Brown #define freq 60
20c4762a1bSJed Brown #define w_s (2 * PETSC_PI * freq)
21c4762a1bSJed Brown
22c4762a1bSJed Brown /* Sizes and indices */
23c4762a1bSJed Brown const PetscInt nbus = 9; /* Number of network buses */
24c4762a1bSJed Brown const PetscInt ngen = 3; /* Number of generators */
25c4762a1bSJed Brown const PetscInt nload = 3; /* Number of loads */
26c4762a1bSJed Brown const PetscInt gbus[3] = {0, 1, 2}; /* Buses at which generators are incident */
27c4762a1bSJed Brown const PetscInt lbus[3] = {4, 5, 7}; /* Buses at which loads are incident */
28c4762a1bSJed Brown
29c4762a1bSJed Brown /* Generator real and reactive powers (found via loadflow) */
30c4762a1bSJed Brown const PetscScalar PG[3] = {0.716786142395021, 1.630000000000000, 0.850000000000000};
31c4762a1bSJed Brown const PetscScalar QG[3] = {0.270702180178785, 0.066120127797275, -0.108402221791588};
32c4762a1bSJed Brown /* Generator constants */
33c4762a1bSJed Brown const PetscScalar H[3] = {23.64, 6.4, 3.01}; /* Inertia constant */
34c4762a1bSJed Brown const PetscScalar Rs[3] = {0.0, 0.0, 0.0}; /* Stator Resistance */
35c4762a1bSJed Brown const PetscScalar Xd[3] = {0.146, 0.8958, 1.3125}; /* d-axis reactance */
36c4762a1bSJed Brown const PetscScalar Xdp[3] = {0.0608, 0.1198, 0.1813}; /* d-axis transient reactance */
37c4762a1bSJed 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 */
38c4762a1bSJed Brown const PetscScalar Xqp[3] = {0.0969, 0.1969, 0.25}; /* q-axis transient reactance */
39c4762a1bSJed Brown const PetscScalar Td0p[3] = {8.96, 6.0, 5.89}; /* d-axis open circuit time constant */
40c4762a1bSJed Brown const PetscScalar Tq0p[3] = {0.31, 0.535, 0.6}; /* q-axis open circuit time constant */
41c4762a1bSJed Brown PetscScalar M[3]; /* M = 2*H/w_s */
42c4762a1bSJed Brown PetscScalar D[3]; /* D = 0.1*M */
43c4762a1bSJed Brown
44c4762a1bSJed Brown PetscScalar TM[3]; /* Mechanical Torque */
45c4762a1bSJed Brown /* Exciter system constants */
46c4762a1bSJed Brown const PetscScalar KA[3] = {20.0, 20.0, 20.0}; /* Voltage regulartor gain constant */
47c4762a1bSJed Brown const PetscScalar TA[3] = {0.2, 0.2, 0.2}; /* Voltage regulator time constant */
48c4762a1bSJed Brown const PetscScalar KE[3] = {1.0, 1.0, 1.0}; /* Exciter gain constant */
49c4762a1bSJed Brown const PetscScalar TE[3] = {0.314, 0.314, 0.314}; /* Exciter time constant */
50c4762a1bSJed Brown const PetscScalar KF[3] = {0.063, 0.063, 0.063}; /* Feedback stabilizer gain constant */
51c4762a1bSJed Brown const PetscScalar TF[3] = {0.35, 0.35, 0.35}; /* Feedback stabilizer time constant */
52c4762a1bSJed Brown const PetscScalar k1[3] = {0.0039, 0.0039, 0.0039};
53c4762a1bSJed Brown const PetscScalar k2[3] = {1.555, 1.555, 1.555}; /* k1 and k2 for calculating the saturation function SE = k1*exp(k2*Efd) */
54c4762a1bSJed Brown
55c4762a1bSJed Brown PetscScalar Vref[3];
56c4762a1bSJed Brown /* Load constants
57c4762a1bSJed Brown We use a composite load model that describes the load and reactive powers at each time instant as follows
58c4762a1bSJed Brown P(t) = \sum\limits_{i=0}^ld_nsegsp \ld_alphap_i*P_D0(\frac{V_m(t)}{V_m0})^\ld_betap_i
59c4762a1bSJed Brown Q(t) = \sum\limits_{i=0}^ld_nsegsq \ld_alphaq_i*Q_D0(\frac{V_m(t)}{V_m0})^\ld_betaq_i
60c4762a1bSJed Brown where
61c4762a1bSJed Brown ld_nsegsp,ld_nsegsq - Number of individual load models for real and reactive power loads
62c4762a1bSJed Brown ld_alphap,ld_alphap - Percentage contribution (weights) or loads
63c4762a1bSJed Brown P_D0 - Real power load
64c4762a1bSJed Brown Q_D0 - Reactive power load
65c4762a1bSJed Brown V_m(t) - Voltage magnitude at time t
66c4762a1bSJed Brown V_m0 - Voltage magnitude at t = 0
67c4762a1bSJed Brown ld_betap, ld_betaq - exponents describing the load model for real and reactive part
68c4762a1bSJed Brown
69c4762a1bSJed Brown Note: All loads have the same characteristic currently.
70c4762a1bSJed Brown */
71c4762a1bSJed Brown const PetscScalar PD0[3] = {1.25, 0.9, 1.0};
72c4762a1bSJed Brown const PetscScalar QD0[3] = {0.5, 0.3, 0.35};
73c4762a1bSJed Brown const PetscInt ld_nsegsp[3] = {3, 3, 3};
74c4762a1bSJed Brown const PetscScalar ld_alphap[3] = {1.0, 0.0, 0.0};
75c4762a1bSJed Brown const PetscScalar ld_betap[3] = {2.0, 1.0, 0.0};
76c4762a1bSJed Brown const PetscInt ld_nsegsq[3] = {3, 3, 3};
77c4762a1bSJed Brown const PetscScalar ld_alphaq[3] = {1.0, 0.0, 0.0};
78c4762a1bSJed Brown const PetscScalar ld_betaq[3] = {2.0, 1.0, 0.0};
79c4762a1bSJed Brown
80c4762a1bSJed Brown typedef struct {
81c4762a1bSJed Brown DM dmgen, dmnet; /* DMs to manage generator and network subsystem */
82c4762a1bSJed Brown DM dmpgrid; /* Composite DM to manage the entire power grid */
83c4762a1bSJed Brown Mat Ybus; /* Network admittance matrix */
84c4762a1bSJed Brown Vec V0; /* Initial voltage vector (Power flow solution) */
85c4762a1bSJed Brown PetscReal tfaulton, tfaultoff; /* Fault on and off times */
86c4762a1bSJed Brown PetscInt faultbus; /* Fault bus */
87c4762a1bSJed Brown PetscScalar Rfault;
88c4762a1bSJed Brown PetscReal t0, tmax;
89c4762a1bSJed Brown PetscInt neqs_gen, neqs_net, neqs_pgrid;
90c4762a1bSJed Brown PetscBool alg_flg;
91c4762a1bSJed Brown PetscReal t;
92c4762a1bSJed Brown IS is_diff; /* indices for differential equations */
93c4762a1bSJed Brown IS is_alg; /* indices for algebraic equations */
94c4762a1bSJed Brown } Userctx;
95c4762a1bSJed Brown
96c4762a1bSJed 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)97d71ae5a4SJacob Faibussowitsch PetscErrorCode dq2ri(PetscScalar Fd, PetscScalar Fq, PetscScalar delta, PetscScalar *Fr, PetscScalar *Fi)
98d71ae5a4SJacob Faibussowitsch {
99c4762a1bSJed Brown PetscFunctionBegin;
100c4762a1bSJed Brown *Fr = Fd * PetscSinScalar(delta) + Fq * PetscCosScalar(delta);
101c4762a1bSJed Brown *Fi = -Fd * PetscCosScalar(delta) + Fq * PetscSinScalar(delta);
1023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
103c4762a1bSJed Brown }
104c4762a1bSJed Brown
105c4762a1bSJed 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)106d71ae5a4SJacob Faibussowitsch PetscErrorCode ri2dq(PetscScalar Fr, PetscScalar Fi, PetscScalar delta, PetscScalar *Fd, PetscScalar *Fq)
107d71ae5a4SJacob Faibussowitsch {
108c4762a1bSJed Brown PetscFunctionBegin;
109c4762a1bSJed Brown *Fd = Fr * PetscSinScalar(delta) - Fi * PetscCosScalar(delta);
110c4762a1bSJed Brown *Fq = Fr * PetscCosScalar(delta) + Fi * PetscSinScalar(delta);
1113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
112c4762a1bSJed Brown }
113c4762a1bSJed Brown
SetInitialGuess(Vec X,Userctx * user)114d71ae5a4SJacob Faibussowitsch PetscErrorCode SetInitialGuess(Vec X, Userctx *user)
115d71ae5a4SJacob Faibussowitsch {
116c4762a1bSJed Brown Vec Xgen, Xnet;
117c4762a1bSJed Brown PetscScalar *xgen, *xnet;
118c4762a1bSJed Brown PetscInt i, idx = 0;
119c4762a1bSJed Brown PetscScalar Vr, Vi, IGr, IGi, Vm, Vm2;
120c4762a1bSJed Brown PetscScalar Eqp, Edp, delta;
121c4762a1bSJed Brown PetscScalar Efd, RF, VR; /* Exciter variables */
122c4762a1bSJed Brown PetscScalar Id, Iq; /* Generator dq axis currents */
123c4762a1bSJed Brown PetscScalar theta, Vd, Vq, SE;
124c4762a1bSJed Brown
125c4762a1bSJed Brown PetscFunctionBegin;
1269371c9d4SSatish Balay M[0] = 2 * H[0] / w_s;
1279371c9d4SSatish Balay M[1] = 2 * H[1] / w_s;
1289371c9d4SSatish Balay M[2] = 2 * H[2] / w_s;
1299371c9d4SSatish Balay D[0] = 0.1 * M[0];
1309371c9d4SSatish Balay D[1] = 0.1 * M[1];
1319371c9d4SSatish Balay D[2] = 0.1 * M[2];
132c4762a1bSJed Brown
1339566063dSJacob Faibussowitsch PetscCall(DMCompositeGetLocalVectors(user->dmpgrid, &Xgen, &Xnet));
134c4762a1bSJed Brown
135c4762a1bSJed Brown /* Network subsystem initialization */
1369566063dSJacob Faibussowitsch PetscCall(VecCopy(user->V0, Xnet));
137c4762a1bSJed Brown
138c4762a1bSJed Brown /* Generator subsystem initialization */
1399566063dSJacob Faibussowitsch PetscCall(VecGetArray(Xgen, &xgen));
1409566063dSJacob Faibussowitsch PetscCall(VecGetArray(Xnet, &xnet));
141c4762a1bSJed Brown
142c4762a1bSJed Brown for (i = 0; i < ngen; i++) {
143c4762a1bSJed Brown Vr = xnet[2 * gbus[i]]; /* Real part of generator terminal voltage */
144c4762a1bSJed Brown Vi = xnet[2 * gbus[i] + 1]; /* Imaginary part of the generator terminal voltage */
1459371c9d4SSatish Balay Vm = PetscSqrtScalar(Vr * Vr + Vi * Vi);
1469371c9d4SSatish Balay Vm2 = Vm * Vm;
147c4762a1bSJed Brown IGr = (Vr * PG[i] + Vi * QG[i]) / Vm2;
148c4762a1bSJed Brown IGi = (Vi * PG[i] - Vr * QG[i]) / Vm2;
149c4762a1bSJed Brown
150c4762a1bSJed Brown delta = PetscAtan2Real(Vi + Xq[i] * IGr, Vr - Xq[i] * IGi); /* Machine angle */
151c4762a1bSJed Brown
152c4762a1bSJed Brown theta = PETSC_PI / 2.0 - delta;
153c4762a1bSJed Brown
154c4762a1bSJed Brown Id = IGr * PetscCosScalar(theta) - IGi * PetscSinScalar(theta); /* d-axis stator current */
155c4762a1bSJed Brown Iq = IGr * PetscSinScalar(theta) + IGi * PetscCosScalar(theta); /* q-axis stator current */
156c4762a1bSJed Brown
157c4762a1bSJed Brown Vd = Vr * PetscCosScalar(theta) - Vi * PetscSinScalar(theta);
158c4762a1bSJed Brown Vq = Vr * PetscSinScalar(theta) + Vi * PetscCosScalar(theta);
159c4762a1bSJed Brown
160c4762a1bSJed Brown Edp = Vd + Rs[i] * Id - Xqp[i] * Iq; /* d-axis transient EMF */
161c4762a1bSJed Brown Eqp = Vq + Rs[i] * Iq + Xdp[i] * Id; /* q-axis transient EMF */
162c4762a1bSJed Brown
163c4762a1bSJed Brown TM[i] = PG[i];
164c4762a1bSJed Brown
165c4762a1bSJed Brown /* The generator variables are ordered as [Eqp,Edp,delta,w,Id,Iq] */
166c4762a1bSJed Brown xgen[idx] = Eqp;
167c4762a1bSJed Brown xgen[idx + 1] = Edp;
168c4762a1bSJed Brown xgen[idx + 2] = delta;
169c4762a1bSJed Brown xgen[idx + 3] = w_s;
170c4762a1bSJed Brown
171c4762a1bSJed Brown idx = idx + 4;
172c4762a1bSJed Brown
173c4762a1bSJed Brown xgen[idx] = Id;
174c4762a1bSJed Brown xgen[idx + 1] = Iq;
175c4762a1bSJed Brown
176c4762a1bSJed Brown idx = idx + 2;
177c4762a1bSJed Brown
178c4762a1bSJed Brown /* Exciter */
179c4762a1bSJed Brown Efd = Eqp + (Xd[i] - Xdp[i]) * Id;
180c4762a1bSJed Brown SE = k1[i] * PetscExpScalar(k2[i] * Efd);
181c4762a1bSJed Brown VR = KE[i] * Efd + SE;
182c4762a1bSJed Brown RF = KF[i] * Efd / TF[i];
183c4762a1bSJed Brown
184c4762a1bSJed Brown xgen[idx] = Efd;
185c4762a1bSJed Brown xgen[idx + 1] = RF;
186c4762a1bSJed Brown xgen[idx + 2] = VR;
187c4762a1bSJed Brown
188c4762a1bSJed Brown Vref[i] = Vm + (VR / KA[i]);
189c4762a1bSJed Brown
190c4762a1bSJed Brown idx = idx + 3;
191c4762a1bSJed Brown }
192c4762a1bSJed Brown
1939566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Xgen, &xgen));
1949566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Xnet, &xnet));
195c4762a1bSJed Brown
1969566063dSJacob Faibussowitsch /* PetscCall(VecView(Xgen,0)); */
1979566063dSJacob Faibussowitsch PetscCall(DMCompositeGather(user->dmpgrid, INSERT_VALUES, X, Xgen, Xnet));
1989566063dSJacob Faibussowitsch PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid, &Xgen, &Xnet));
1993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
200c4762a1bSJed Brown }
201c4762a1bSJed Brown
202c4762a1bSJed Brown /* Computes F = [f(x,y);g(x,y)] */
ResidualFunction(SNES snes,Vec X,Vec F,Userctx * user)203d71ae5a4SJacob Faibussowitsch PetscErrorCode ResidualFunction(SNES snes, Vec X, Vec F, Userctx *user)
204d71ae5a4SJacob Faibussowitsch {
205c4762a1bSJed Brown Vec Xgen, Xnet, Fgen, Fnet;
206c4762a1bSJed Brown PetscScalar *xgen, *xnet, *fgen, *fnet;
207c4762a1bSJed Brown PetscInt i, idx = 0;
208c4762a1bSJed Brown PetscScalar Vr, Vi, Vm, Vm2;
209c4762a1bSJed Brown PetscScalar Eqp, Edp, delta, w; /* Generator variables */
210c4762a1bSJed Brown PetscScalar Efd, RF, VR; /* Exciter variables */
211c4762a1bSJed Brown PetscScalar Id, Iq; /* Generator dq axis currents */
212c4762a1bSJed Brown PetscScalar Vd, Vq, SE;
213c4762a1bSJed Brown PetscScalar IGr, IGi, IDr, IDi;
214c4762a1bSJed Brown PetscScalar Zdq_inv[4], det;
215c4762a1bSJed Brown PetscScalar PD, QD, Vm0, *v0;
216c4762a1bSJed Brown PetscInt k;
217c4762a1bSJed Brown
218c4762a1bSJed Brown PetscFunctionBegin;
2199566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F));
2209566063dSJacob Faibussowitsch PetscCall(DMCompositeGetLocalVectors(user->dmpgrid, &Xgen, &Xnet));
2219566063dSJacob Faibussowitsch PetscCall(DMCompositeGetLocalVectors(user->dmpgrid, &Fgen, &Fnet));
2229566063dSJacob Faibussowitsch PetscCall(DMCompositeScatter(user->dmpgrid, X, Xgen, Xnet));
2239566063dSJacob Faibussowitsch PetscCall(DMCompositeScatter(user->dmpgrid, F, Fgen, Fnet));
224c4762a1bSJed Brown
225c4762a1bSJed Brown /* Network current balance residual IG + Y*V + IL = 0. Only YV is added here.
226c4762a1bSJed Brown The generator current injection, IG, and load current injection, ID are added later
227c4762a1bSJed Brown */
228c4762a1bSJed Brown /* Note that the values in Ybus are stored assuming the imaginary current balance
229c4762a1bSJed Brown equation is ordered first followed by real current balance equation for each bus.
230c4762a1bSJed Brown Thus imaginary current contribution goes in location 2*i, and
231c4762a1bSJed Brown real current contribution in 2*i+1
232c4762a1bSJed Brown */
2339566063dSJacob Faibussowitsch PetscCall(MatMult(user->Ybus, Xnet, Fnet));
234c4762a1bSJed Brown
2359566063dSJacob Faibussowitsch PetscCall(VecGetArray(Xgen, &xgen));
2369566063dSJacob Faibussowitsch PetscCall(VecGetArray(Xnet, &xnet));
2379566063dSJacob Faibussowitsch PetscCall(VecGetArray(Fgen, &fgen));
2389566063dSJacob Faibussowitsch PetscCall(VecGetArray(Fnet, &fnet));
239c4762a1bSJed Brown
240c4762a1bSJed Brown /* Generator subsystem */
241c4762a1bSJed Brown for (i = 0; i < ngen; i++) {
242c4762a1bSJed Brown Eqp = xgen[idx];
243c4762a1bSJed Brown Edp = xgen[idx + 1];
244c4762a1bSJed Brown delta = xgen[idx + 2];
245c4762a1bSJed Brown w = xgen[idx + 3];
246c4762a1bSJed Brown Id = xgen[idx + 4];
247c4762a1bSJed Brown Iq = xgen[idx + 5];
248c4762a1bSJed Brown Efd = xgen[idx + 6];
249c4762a1bSJed Brown RF = xgen[idx + 7];
250c4762a1bSJed Brown VR = xgen[idx + 8];
251c4762a1bSJed Brown
252c4762a1bSJed Brown /* Generator differential equations */
253c4762a1bSJed Brown fgen[idx] = (Eqp + (Xd[i] - Xdp[i]) * Id - Efd) / Td0p[i];
254c4762a1bSJed Brown fgen[idx + 1] = (Edp - (Xq[i] - Xqp[i]) * Iq) / Tq0p[i];
255c4762a1bSJed Brown fgen[idx + 2] = -w + w_s;
256c4762a1bSJed Brown fgen[idx + 3] = (-TM[i] + Edp * Id + Eqp * Iq + (Xqp[i] - Xdp[i]) * Id * Iq + D[i] * (w - w_s)) / M[i];
257c4762a1bSJed Brown
258c4762a1bSJed Brown Vr = xnet[2 * gbus[i]]; /* Real part of generator terminal voltage */
259c4762a1bSJed Brown Vi = xnet[2 * gbus[i] + 1]; /* Imaginary part of the generator terminal voltage */
260c4762a1bSJed Brown
2619566063dSJacob Faibussowitsch PetscCall(ri2dq(Vr, Vi, delta, &Vd, &Vq));
262c4762a1bSJed Brown /* Algebraic equations for stator currents */
263c4762a1bSJed Brown
264c4762a1bSJed Brown det = Rs[i] * Rs[i] + Xdp[i] * Xqp[i];
265c4762a1bSJed Brown
266c4762a1bSJed Brown Zdq_inv[0] = Rs[i] / det;
267c4762a1bSJed Brown Zdq_inv[1] = Xqp[i] / det;
268c4762a1bSJed Brown Zdq_inv[2] = -Xdp[i] / det;
269c4762a1bSJed Brown Zdq_inv[3] = Rs[i] / det;
270c4762a1bSJed Brown
271c4762a1bSJed Brown fgen[idx + 4] = Zdq_inv[0] * (-Edp + Vd) + Zdq_inv[1] * (-Eqp + Vq) + Id;
272c4762a1bSJed Brown fgen[idx + 5] = Zdq_inv[2] * (-Edp + Vd) + Zdq_inv[3] * (-Eqp + Vq) + Iq;
273c4762a1bSJed Brown
274c4762a1bSJed Brown /* Add generator current injection to network */
2759566063dSJacob Faibussowitsch PetscCall(dq2ri(Id, Iq, delta, &IGr, &IGi));
276c4762a1bSJed Brown
277c4762a1bSJed Brown fnet[2 * gbus[i]] -= IGi;
278c4762a1bSJed Brown fnet[2 * gbus[i] + 1] -= IGr;
279c4762a1bSJed Brown
280c4762a1bSJed Brown Vm = PetscSqrtScalar(Vd * Vd + Vq * Vq);
281c4762a1bSJed Brown
282c4762a1bSJed Brown SE = k1[i] * PetscExpScalar(k2[i] * Efd);
283c4762a1bSJed Brown
284c4762a1bSJed Brown /* Exciter differential equations */
285c4762a1bSJed Brown fgen[idx + 6] = (KE[i] * Efd + SE - VR) / TE[i];
286c4762a1bSJed Brown fgen[idx + 7] = (RF - KF[i] * Efd / TF[i]) / TF[i];
287c4762a1bSJed Brown fgen[idx + 8] = (VR - KA[i] * RF + KA[i] * KF[i] * Efd / TF[i] - KA[i] * (Vref[i] - Vm)) / TA[i];
288c4762a1bSJed Brown
289c4762a1bSJed Brown idx = idx + 9;
290c4762a1bSJed Brown }
291c4762a1bSJed Brown
2929566063dSJacob Faibussowitsch PetscCall(VecGetArray(user->V0, &v0));
293c4762a1bSJed Brown for (i = 0; i < nload; i++) {
294c4762a1bSJed Brown Vr = xnet[2 * lbus[i]]; /* Real part of load bus voltage */
295c4762a1bSJed Brown Vi = xnet[2 * lbus[i] + 1]; /* Imaginary part of the load bus voltage */
2969371c9d4SSatish Balay Vm = PetscSqrtScalar(Vr * Vr + Vi * Vi);
2979371c9d4SSatish Balay Vm2 = Vm * Vm;
298c4762a1bSJed Brown Vm0 = PetscSqrtScalar(v0[2 * lbus[i]] * v0[2 * lbus[i]] + v0[2 * lbus[i] + 1] * v0[2 * lbus[i] + 1]);
299c4762a1bSJed Brown PD = QD = 0.0;
30057508eceSPierre Jolivet for (k = 0; k < ld_nsegsp[i]; k++) PD += ld_alphap[k] * PD0[i] * PetscPowScalar(Vm / Vm0, ld_betap[k]);
30157508eceSPierre Jolivet for (k = 0; k < ld_nsegsq[i]; k++) QD += ld_alphaq[k] * QD0[i] * PetscPowScalar(Vm / Vm0, ld_betaq[k]);
302c4762a1bSJed Brown
303c4762a1bSJed Brown /* Load currents */
304c4762a1bSJed Brown IDr = (PD * Vr + QD * Vi) / Vm2;
305c4762a1bSJed Brown IDi = (-QD * Vr + PD * Vi) / Vm2;
306c4762a1bSJed Brown
307c4762a1bSJed Brown fnet[2 * lbus[i]] += IDi;
308c4762a1bSJed Brown fnet[2 * lbus[i] + 1] += IDr;
309c4762a1bSJed Brown }
3109566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(user->V0, &v0));
311c4762a1bSJed Brown
3129566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Xgen, &xgen));
3139566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Xnet, &xnet));
3149566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Fgen, &fgen));
3159566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Fnet, &fnet));
316c4762a1bSJed Brown
3179566063dSJacob Faibussowitsch PetscCall(DMCompositeGather(user->dmpgrid, INSERT_VALUES, F, Fgen, Fnet));
3189566063dSJacob Faibussowitsch PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid, &Xgen, &Xnet));
3199566063dSJacob Faibussowitsch PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid, &Fgen, &Fnet));
3203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
321c4762a1bSJed Brown }
322c4762a1bSJed Brown
323c4762a1bSJed Brown /* \dot{x} - f(x,y)
324c4762a1bSJed Brown g(x,y) = 0
325c4762a1bSJed Brown */
IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,Userctx * user)326d71ae5a4SJacob Faibussowitsch PetscErrorCode IFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, Userctx *user)
327d71ae5a4SJacob Faibussowitsch {
328c4762a1bSJed Brown SNES snes;
329c4762a1bSJed Brown PetscScalar *f;
330c4762a1bSJed Brown const PetscScalar *xdot;
331c4762a1bSJed Brown PetscInt i;
332c4762a1bSJed Brown
333c4762a1bSJed Brown PetscFunctionBegin;
334c4762a1bSJed Brown user->t = t;
335c4762a1bSJed Brown
3369566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes));
3379566063dSJacob Faibussowitsch PetscCall(ResidualFunction(snes, X, F, user));
3389566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f));
3399566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Xdot, &xdot));
340c4762a1bSJed Brown for (i = 0; i < ngen; i++) {
341c4762a1bSJed Brown f[9 * i] += xdot[9 * i];
342c4762a1bSJed Brown f[9 * i + 1] += xdot[9 * i + 1];
343c4762a1bSJed Brown f[9 * i + 2] += xdot[9 * i + 2];
344c4762a1bSJed Brown f[9 * i + 3] += xdot[9 * i + 3];
345c4762a1bSJed Brown f[9 * i + 6] += xdot[9 * i + 6];
346c4762a1bSJed Brown f[9 * i + 7] += xdot[9 * i + 7];
347c4762a1bSJed Brown f[9 * i + 8] += xdot[9 * i + 8];
348c4762a1bSJed Brown }
3499566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f));
3509566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Xdot, &xdot));
3513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
352c4762a1bSJed Brown }
353c4762a1bSJed Brown
354c4762a1bSJed Brown /* This function is used for solving the algebraic system only during fault on and
355c4762a1bSJed Brown off times. It computes the entire F and then zeros out the part corresponding to
356c4762a1bSJed Brown differential equations
357c4762a1bSJed Brown F = [0;g(y)];
358c4762a1bSJed Brown */
AlgFunction(SNES snes,Vec X,Vec F,PetscCtx ctx)359*2a8381b2SBarry Smith PetscErrorCode AlgFunction(SNES snes, Vec X, Vec F, PetscCtx ctx)
360d71ae5a4SJacob Faibussowitsch {
361c4762a1bSJed Brown Userctx *user = (Userctx *)ctx;
362c4762a1bSJed Brown PetscScalar *f;
363c4762a1bSJed Brown PetscInt i;
364c4762a1bSJed Brown
365c4762a1bSJed Brown PetscFunctionBegin;
3669566063dSJacob Faibussowitsch PetscCall(ResidualFunction(snes, X, F, user));
3679566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f));
368c4762a1bSJed Brown for (i = 0; i < ngen; i++) {
369c4762a1bSJed Brown f[9 * i] = 0;
370c4762a1bSJed Brown f[9 * i + 1] = 0;
371c4762a1bSJed Brown f[9 * i + 2] = 0;
372c4762a1bSJed Brown f[9 * i + 3] = 0;
373c4762a1bSJed Brown f[9 * i + 6] = 0;
374c4762a1bSJed Brown f[9 * i + 7] = 0;
375c4762a1bSJed Brown f[9 * i + 8] = 0;
376c4762a1bSJed Brown }
3779566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f));
3783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
379c4762a1bSJed Brown }
380c4762a1bSJed Brown
PreallocateJacobian(Mat J,Userctx * user)381d71ae5a4SJacob Faibussowitsch PetscErrorCode PreallocateJacobian(Mat J, Userctx *user)
382d71ae5a4SJacob Faibussowitsch {
383c4762a1bSJed Brown PetscInt *d_nnz;
384c4762a1bSJed Brown PetscInt i, idx = 0, start = 0;
385c4762a1bSJed Brown PetscInt ncols;
386c4762a1bSJed Brown
387c4762a1bSJed Brown PetscFunctionBegin;
3889566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->neqs_pgrid, &d_nnz));
389c4762a1bSJed Brown for (i = 0; i < user->neqs_pgrid; i++) d_nnz[i] = 0;
390c4762a1bSJed Brown /* Generator subsystem */
391c4762a1bSJed Brown for (i = 0; i < ngen; i++) {
392c4762a1bSJed Brown d_nnz[idx] += 3;
393c4762a1bSJed Brown d_nnz[idx + 1] += 2;
394c4762a1bSJed Brown d_nnz[idx + 2] += 2;
395c4762a1bSJed Brown d_nnz[idx + 3] += 5;
396c4762a1bSJed Brown d_nnz[idx + 4] += 6;
397c4762a1bSJed Brown d_nnz[idx + 5] += 6;
398c4762a1bSJed Brown
399c4762a1bSJed Brown d_nnz[user->neqs_gen + 2 * gbus[i]] += 3;
400c4762a1bSJed Brown d_nnz[user->neqs_gen + 2 * gbus[i] + 1] += 3;
401c4762a1bSJed Brown
402c4762a1bSJed Brown d_nnz[idx + 6] += 2;
403c4762a1bSJed Brown d_nnz[idx + 7] += 2;
404c4762a1bSJed Brown d_nnz[idx + 8] += 5;
405c4762a1bSJed Brown
406c4762a1bSJed Brown idx = idx + 9;
407c4762a1bSJed Brown }
408c4762a1bSJed Brown
409c4762a1bSJed Brown start = user->neqs_gen;
410c4762a1bSJed Brown
411c4762a1bSJed Brown for (i = 0; i < nbus; i++) {
4129566063dSJacob Faibussowitsch PetscCall(MatGetRow(user->Ybus, 2 * i, &ncols, NULL, NULL));
413c4762a1bSJed Brown d_nnz[start + 2 * i] += ncols;
414c4762a1bSJed Brown d_nnz[start + 2 * i + 1] += ncols;
4159566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(user->Ybus, 2 * i, &ncols, NULL, NULL));
416c4762a1bSJed Brown }
417c4762a1bSJed Brown
4189566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(J, 0, d_nnz));
419c4762a1bSJed Brown
4209566063dSJacob Faibussowitsch PetscCall(PetscFree(d_nnz));
4213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
422c4762a1bSJed Brown }
423c4762a1bSJed Brown
424c4762a1bSJed Brown /*
425c4762a1bSJed Brown J = [-df_dx, -df_dy
426c4762a1bSJed Brown dg_dx, dg_dy]
427c4762a1bSJed Brown */
ResidualJacobian(SNES snes,Vec X,Mat J,Mat B,PetscCtx ctx)428*2a8381b2SBarry Smith PetscErrorCode ResidualJacobian(SNES snes, Vec X, Mat J, Mat B, PetscCtx ctx)
429d71ae5a4SJacob Faibussowitsch {
430c4762a1bSJed Brown Userctx *user = (Userctx *)ctx;
431c4762a1bSJed Brown Vec Xgen, Xnet;
432c4762a1bSJed Brown PetscScalar *xgen, *xnet;
433c4762a1bSJed Brown PetscInt i, idx = 0;
434c4762a1bSJed Brown PetscScalar Vr, Vi, Vm, Vm2;
435c4762a1bSJed Brown PetscScalar Eqp, Edp, delta; /* Generator variables */
436c4762a1bSJed Brown PetscScalar Efd; /* Exciter variables */
437c4762a1bSJed Brown PetscScalar Id, Iq; /* Generator dq axis currents */
438c4762a1bSJed Brown PetscScalar Vd, Vq;
439c4762a1bSJed Brown PetscScalar val[10];
440c4762a1bSJed Brown PetscInt row[2], col[10];
441c4762a1bSJed Brown PetscInt net_start = user->neqs_gen;
442c4762a1bSJed Brown PetscScalar Zdq_inv[4], det;
443c4762a1bSJed Brown PetscScalar dVd_dVr, dVd_dVi, dVq_dVr, dVq_dVi, dVd_ddelta, dVq_ddelta;
444c4762a1bSJed Brown PetscScalar dIGr_ddelta, dIGi_ddelta, dIGr_dId, dIGr_dIq, dIGi_dId, dIGi_dIq;
445c4762a1bSJed Brown PetscScalar dSE_dEfd;
446c4762a1bSJed Brown PetscScalar dVm_dVd, dVm_dVq, dVm_dVr, dVm_dVi;
447c4762a1bSJed Brown PetscInt ncols;
448c4762a1bSJed Brown const PetscInt *cols;
449c4762a1bSJed Brown const PetscScalar *yvals;
450c4762a1bSJed Brown PetscInt k;
451c4762a1bSJed Brown PetscScalar PD, QD, Vm0, *v0, Vm4;
452c4762a1bSJed Brown PetscScalar dPD_dVr, dPD_dVi, dQD_dVr, dQD_dVi;
453c4762a1bSJed Brown PetscScalar dIDr_dVr, dIDr_dVi, dIDi_dVr, dIDi_dVi;
454c4762a1bSJed Brown
455c4762a1bSJed Brown PetscFunctionBegin;
4569566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(B));
4579566063dSJacob Faibussowitsch PetscCall(DMCompositeGetLocalVectors(user->dmpgrid, &Xgen, &Xnet));
4589566063dSJacob Faibussowitsch PetscCall(DMCompositeScatter(user->dmpgrid, X, Xgen, Xnet));
459c4762a1bSJed Brown
4609566063dSJacob Faibussowitsch PetscCall(VecGetArray(Xgen, &xgen));
4619566063dSJacob Faibussowitsch PetscCall(VecGetArray(Xnet, &xnet));
462c4762a1bSJed Brown
463c4762a1bSJed Brown /* Generator subsystem */
464c4762a1bSJed Brown for (i = 0; i < ngen; i++) {
465c4762a1bSJed Brown Eqp = xgen[idx];
466c4762a1bSJed Brown Edp = xgen[idx + 1];
467c4762a1bSJed Brown delta = xgen[idx + 2];
468c4762a1bSJed Brown Id = xgen[idx + 4];
469c4762a1bSJed Brown Iq = xgen[idx + 5];
470c4762a1bSJed Brown Efd = xgen[idx + 6];
471c4762a1bSJed Brown
472c4762a1bSJed Brown /* fgen[idx] = (Eqp + (Xd[i] - Xdp[i])*Id - Efd)/Td0p[i]; */
473c4762a1bSJed Brown row[0] = idx;
4749371c9d4SSatish Balay col[0] = idx;
4759371c9d4SSatish Balay col[1] = idx + 4;
4769371c9d4SSatish Balay col[2] = idx + 6;
4779371c9d4SSatish Balay val[0] = 1 / Td0p[i];
4789371c9d4SSatish Balay val[1] = (Xd[i] - Xdp[i]) / Td0p[i];
4799371c9d4SSatish Balay val[2] = -1 / Td0p[i];
480c4762a1bSJed Brown
4819566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, row, 3, col, val, INSERT_VALUES));
482c4762a1bSJed Brown
483c4762a1bSJed Brown /* fgen[idx+1] = (Edp - (Xq[i] - Xqp[i])*Iq)/Tq0p[i]; */
484c4762a1bSJed Brown row[0] = idx + 1;
4859371c9d4SSatish Balay col[0] = idx + 1;
4869371c9d4SSatish Balay col[1] = idx + 5;
4879371c9d4SSatish Balay val[0] = 1 / Tq0p[i];
4889371c9d4SSatish Balay val[1] = -(Xq[i] - Xqp[i]) / Tq0p[i];
4899566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, row, 2, col, val, INSERT_VALUES));
490c4762a1bSJed Brown
491c4762a1bSJed Brown /* fgen[idx+2] = - w + w_s; */
492c4762a1bSJed Brown row[0] = idx + 2;
4939371c9d4SSatish Balay col[0] = idx + 2;
4949371c9d4SSatish Balay col[1] = idx + 3;
4959371c9d4SSatish Balay val[0] = 0;
4969371c9d4SSatish Balay val[1] = -1;
4979566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, row, 2, col, val, INSERT_VALUES));
498c4762a1bSJed Brown
499c4762a1bSJed Brown /* fgen[idx+3] = (-TM[i] + Edp*Id + Eqp*Iq + (Xqp[i] - Xdp[i])*Id*Iq + D[i]*(w - w_s))/M[i]; */
500c4762a1bSJed Brown row[0] = idx + 3;
5019371c9d4SSatish Balay col[0] = idx;
5029371c9d4SSatish Balay col[1] = idx + 1;
5039371c9d4SSatish Balay col[2] = idx + 3;
5049371c9d4SSatish Balay col[3] = idx + 4;
5059371c9d4SSatish Balay col[4] = idx + 5;
5069371c9d4SSatish Balay val[0] = Iq / M[i];
5079371c9d4SSatish Balay val[1] = Id / M[i];
5089371c9d4SSatish Balay val[2] = D[i] / M[i];
5099371c9d4SSatish Balay val[3] = (Edp + (Xqp[i] - Xdp[i]) * Iq) / M[i];
5109371c9d4SSatish Balay val[4] = (Eqp + (Xqp[i] - Xdp[i]) * Id) / M[i];
5119566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, row, 5, col, val, INSERT_VALUES));
512c4762a1bSJed Brown
513c4762a1bSJed Brown Vr = xnet[2 * gbus[i]]; /* Real part of generator terminal voltage */
514c4762a1bSJed Brown Vi = xnet[2 * gbus[i] + 1]; /* Imaginary part of the generator terminal voltage */
5159566063dSJacob Faibussowitsch PetscCall(ri2dq(Vr, Vi, delta, &Vd, &Vq));
516c4762a1bSJed Brown
517c4762a1bSJed Brown det = Rs[i] * Rs[i] + Xdp[i] * Xqp[i];
518c4762a1bSJed Brown
519c4762a1bSJed Brown Zdq_inv[0] = Rs[i] / det;
520c4762a1bSJed Brown Zdq_inv[1] = Xqp[i] / det;
521c4762a1bSJed Brown Zdq_inv[2] = -Xdp[i] / det;
522c4762a1bSJed Brown Zdq_inv[3] = Rs[i] / det;
523c4762a1bSJed Brown
5249371c9d4SSatish Balay dVd_dVr = PetscSinScalar(delta);
5259371c9d4SSatish Balay dVd_dVi = -PetscCosScalar(delta);
5269371c9d4SSatish Balay dVq_dVr = PetscCosScalar(delta);
5279371c9d4SSatish Balay dVq_dVi = PetscSinScalar(delta);
528c4762a1bSJed Brown dVd_ddelta = Vr * PetscCosScalar(delta) + Vi * PetscSinScalar(delta);
529c4762a1bSJed Brown dVq_ddelta = -Vr * PetscSinScalar(delta) + Vi * PetscCosScalar(delta);
530c4762a1bSJed Brown
531c4762a1bSJed Brown /* fgen[idx+4] = Zdq_inv[0]*(-Edp + Vd) + Zdq_inv[1]*(-Eqp + Vq) + Id; */
532c4762a1bSJed Brown row[0] = idx + 4;
5339371c9d4SSatish Balay col[0] = idx;
5349371c9d4SSatish Balay col[1] = idx + 1;
5359371c9d4SSatish Balay col[2] = idx + 2;
5369371c9d4SSatish Balay val[0] = -Zdq_inv[1];
5379371c9d4SSatish Balay val[1] = -Zdq_inv[0];
5389371c9d4SSatish Balay val[2] = Zdq_inv[0] * dVd_ddelta + Zdq_inv[1] * dVq_ddelta;
5399371c9d4SSatish Balay col[3] = idx + 4;
5409371c9d4SSatish Balay col[4] = net_start + 2 * gbus[i];
5419371c9d4SSatish Balay col[5] = net_start + 2 * gbus[i] + 1;
5429371c9d4SSatish Balay val[3] = 1;
5439371c9d4SSatish Balay val[4] = Zdq_inv[0] * dVd_dVr + Zdq_inv[1] * dVq_dVr;
5449371c9d4SSatish Balay val[5] = Zdq_inv[0] * dVd_dVi + Zdq_inv[1] * dVq_dVi;
5459566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, row, 6, col, val, INSERT_VALUES));
546c4762a1bSJed Brown
547c4762a1bSJed Brown /* fgen[idx+5] = Zdq_inv[2]*(-Edp + Vd) + Zdq_inv[3]*(-Eqp + Vq) + Iq; */
548c4762a1bSJed Brown row[0] = idx + 5;
5499371c9d4SSatish Balay col[0] = idx;
5509371c9d4SSatish Balay col[1] = idx + 1;
5519371c9d4SSatish Balay col[2] = idx + 2;
5529371c9d4SSatish Balay val[0] = -Zdq_inv[3];
5539371c9d4SSatish Balay val[1] = -Zdq_inv[2];
5549371c9d4SSatish Balay val[2] = Zdq_inv[2] * dVd_ddelta + Zdq_inv[3] * dVq_ddelta;
5559371c9d4SSatish Balay col[3] = idx + 5;
5569371c9d4SSatish Balay col[4] = net_start + 2 * gbus[i];
5579371c9d4SSatish Balay col[5] = net_start + 2 * gbus[i] + 1;
5589371c9d4SSatish Balay val[3] = 1;
5599371c9d4SSatish Balay val[4] = Zdq_inv[2] * dVd_dVr + Zdq_inv[3] * dVq_dVr;
5609371c9d4SSatish Balay val[5] = Zdq_inv[2] * dVd_dVi + Zdq_inv[3] * dVq_dVi;
5619566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, row, 6, col, val, INSERT_VALUES));
562c4762a1bSJed Brown
563c4762a1bSJed Brown dIGr_ddelta = Id * PetscCosScalar(delta) - Iq * PetscSinScalar(delta);
564c4762a1bSJed Brown dIGi_ddelta = Id * PetscSinScalar(delta) + Iq * PetscCosScalar(delta);
5659371c9d4SSatish Balay dIGr_dId = PetscSinScalar(delta);
5669371c9d4SSatish Balay dIGr_dIq = PetscCosScalar(delta);
5679371c9d4SSatish Balay dIGi_dId = -PetscCosScalar(delta);
5689371c9d4SSatish Balay dIGi_dIq = PetscSinScalar(delta);
569c4762a1bSJed Brown
570c4762a1bSJed Brown /* fnet[2*gbus[i]] -= IGi; */
571c4762a1bSJed Brown row[0] = net_start + 2 * gbus[i];
5729371c9d4SSatish Balay col[0] = idx + 2;
5739371c9d4SSatish Balay col[1] = idx + 4;
5749371c9d4SSatish Balay col[2] = idx + 5;
5759371c9d4SSatish Balay val[0] = -dIGi_ddelta;
5769371c9d4SSatish Balay val[1] = -dIGi_dId;
5779371c9d4SSatish Balay val[2] = -dIGi_dIq;
5789566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, row, 3, col, val, INSERT_VALUES));
579c4762a1bSJed Brown
580c4762a1bSJed Brown /* fnet[2*gbus[i]+1] -= IGr; */
581c4762a1bSJed Brown row[0] = net_start + 2 * gbus[i] + 1;
5829371c9d4SSatish Balay col[0] = idx + 2;
5839371c9d4SSatish Balay col[1] = idx + 4;
5849371c9d4SSatish Balay col[2] = idx + 5;
5859371c9d4SSatish Balay val[0] = -dIGr_ddelta;
5869371c9d4SSatish Balay val[1] = -dIGr_dId;
5879371c9d4SSatish Balay val[2] = -dIGr_dIq;
5889566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, row, 3, col, val, INSERT_VALUES));
589c4762a1bSJed Brown
590c4762a1bSJed Brown Vm = PetscSqrtScalar(Vd * Vd + Vq * Vq);
591c4762a1bSJed Brown
592c4762a1bSJed Brown /* fgen[idx+6] = (KE[i]*Efd + SE - VR)/TE[i]; */
593c4762a1bSJed Brown /* SE = k1[i]*PetscExpScalar(k2[i]*Efd); */
594c4762a1bSJed Brown dSE_dEfd = k1[i] * k2[i] * PetscExpScalar(k2[i] * Efd);
595c4762a1bSJed Brown
596c4762a1bSJed Brown row[0] = idx + 6;
5979371c9d4SSatish Balay col[0] = idx + 6;
5989371c9d4SSatish Balay col[1] = idx + 8;
5999371c9d4SSatish Balay val[0] = (KE[i] + dSE_dEfd) / TE[i];
6009371c9d4SSatish Balay val[1] = -1 / TE[i];
6019566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, row, 2, col, val, INSERT_VALUES));
602c4762a1bSJed Brown
603c4762a1bSJed Brown /* Exciter differential equations */
604c4762a1bSJed Brown
605c4762a1bSJed Brown /* fgen[idx+7] = (RF - KF[i]*Efd/TF[i])/TF[i]; */
606c4762a1bSJed Brown row[0] = idx + 7;
6079371c9d4SSatish Balay col[0] = idx + 6;
6089371c9d4SSatish Balay col[1] = idx + 7;
6099371c9d4SSatish Balay val[0] = (-KF[i] / TF[i]) / TF[i];
6109371c9d4SSatish Balay val[1] = 1 / TF[i];
6119566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, row, 2, col, val, INSERT_VALUES));
612c4762a1bSJed Brown
613c4762a1bSJed Brown /* fgen[idx+8] = (VR - KA[i]*RF + KA[i]*KF[i]*Efd/TF[i] - KA[i]*(Vref[i] - Vm))/TA[i]; */
614c4762a1bSJed Brown /* Vm = (Vd^2 + Vq^2)^0.5; */
6159371c9d4SSatish Balay dVm_dVd = Vd / Vm;
6169371c9d4SSatish Balay dVm_dVq = Vq / Vm;
617c4762a1bSJed Brown dVm_dVr = dVm_dVd * dVd_dVr + dVm_dVq * dVq_dVr;
618c4762a1bSJed Brown dVm_dVi = dVm_dVd * dVd_dVi + dVm_dVq * dVq_dVi;
619c4762a1bSJed Brown row[0] = idx + 8;
6209371c9d4SSatish Balay col[0] = idx + 6;
6219371c9d4SSatish Balay col[1] = idx + 7;
6229371c9d4SSatish Balay col[2] = idx + 8;
6239371c9d4SSatish Balay val[0] = (KA[i] * KF[i] / TF[i]) / TA[i];
6249371c9d4SSatish Balay val[1] = -KA[i] / TA[i];
6259371c9d4SSatish Balay val[2] = 1 / TA[i];
6269371c9d4SSatish Balay col[3] = net_start + 2 * gbus[i];
6279371c9d4SSatish Balay col[4] = net_start + 2 * gbus[i] + 1;
6289371c9d4SSatish Balay val[3] = KA[i] * dVm_dVr / TA[i];
6299371c9d4SSatish Balay val[4] = KA[i] * dVm_dVi / TA[i];
6309566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, row, 5, col, val, INSERT_VALUES));
631c4762a1bSJed Brown idx = idx + 9;
632c4762a1bSJed Brown }
633c4762a1bSJed Brown
634c4762a1bSJed Brown for (i = 0; i < nbus; i++) {
6359566063dSJacob Faibussowitsch PetscCall(MatGetRow(user->Ybus, 2 * i, &ncols, &cols, &yvals));
636c4762a1bSJed Brown row[0] = net_start + 2 * i;
637c4762a1bSJed Brown for (k = 0; k < ncols; k++) {
638c4762a1bSJed Brown col[k] = net_start + cols[k];
639c4762a1bSJed Brown val[k] = yvals[k];
640c4762a1bSJed Brown }
6419566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, row, ncols, col, val, INSERT_VALUES));
6429566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(user->Ybus, 2 * i, &ncols, &cols, &yvals));
643c4762a1bSJed Brown
6449566063dSJacob Faibussowitsch PetscCall(MatGetRow(user->Ybus, 2 * i + 1, &ncols, &cols, &yvals));
645c4762a1bSJed Brown row[0] = net_start + 2 * i + 1;
646c4762a1bSJed Brown for (k = 0; k < ncols; k++) {
647c4762a1bSJed Brown col[k] = net_start + cols[k];
648c4762a1bSJed Brown val[k] = yvals[k];
649c4762a1bSJed Brown }
6509566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, row, ncols, col, val, INSERT_VALUES));
6519566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(user->Ybus, 2 * i + 1, &ncols, &cols, &yvals));
652c4762a1bSJed Brown }
653c4762a1bSJed Brown
6549566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J, MAT_FLUSH_ASSEMBLY));
6559566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J, MAT_FLUSH_ASSEMBLY));
656c4762a1bSJed Brown
6579566063dSJacob Faibussowitsch PetscCall(VecGetArray(user->V0, &v0));
658c4762a1bSJed Brown for (i = 0; i < nload; i++) {
659c4762a1bSJed Brown Vr = xnet[2 * lbus[i]]; /* Real part of load bus voltage */
660c4762a1bSJed Brown Vi = xnet[2 * lbus[i] + 1]; /* Imaginary part of the load bus voltage */
6619371c9d4SSatish Balay Vm = PetscSqrtScalar(Vr * Vr + Vi * Vi);
6629371c9d4SSatish Balay Vm2 = Vm * Vm;
6639371c9d4SSatish Balay Vm4 = Vm2 * Vm2;
664c4762a1bSJed Brown Vm0 = PetscSqrtScalar(v0[2 * lbus[i]] * v0[2 * lbus[i]] + v0[2 * lbus[i] + 1] * v0[2 * lbus[i] + 1]);
665c4762a1bSJed Brown PD = QD = 0.0;
666c4762a1bSJed Brown dPD_dVr = dPD_dVi = dQD_dVr = dQD_dVi = 0.0;
667c4762a1bSJed Brown for (k = 0; k < ld_nsegsp[i]; k++) {
66857508eceSPierre Jolivet PD += ld_alphap[k] * PD0[i] * PetscPowScalar(Vm / Vm0, ld_betap[k]);
66957508eceSPierre Jolivet dPD_dVr += ld_alphap[k] * ld_betap[k] * PD0[i] * PetscPowScalar(1 / Vm0, ld_betap[k]) * Vr * PetscPowScalar(Vm, ld_betap[k] - 2);
67057508eceSPierre Jolivet dPD_dVi += ld_alphap[k] * ld_betap[k] * PD0[i] * PetscPowScalar(1 / Vm0, ld_betap[k]) * Vi * PetscPowScalar(Vm, ld_betap[k] - 2);
671c4762a1bSJed Brown }
672c4762a1bSJed Brown for (k = 0; k < ld_nsegsq[i]; k++) {
67357508eceSPierre Jolivet QD += ld_alphaq[k] * QD0[i] * PetscPowScalar(Vm / Vm0, ld_betaq[k]);
67457508eceSPierre Jolivet dQD_dVr += ld_alphaq[k] * ld_betaq[k] * QD0[i] * PetscPowScalar(1 / Vm0, ld_betaq[k]) * Vr * PetscPowScalar(Vm, ld_betaq[k] - 2);
67557508eceSPierre Jolivet dQD_dVi += ld_alphaq[k] * ld_betaq[k] * QD0[i] * PetscPowScalar(1 / Vm0, ld_betaq[k]) * Vi * PetscPowScalar(Vm, ld_betaq[k] - 2);
676c4762a1bSJed Brown }
677c4762a1bSJed Brown
678c4762a1bSJed Brown /* IDr = (PD*Vr + QD*Vi)/Vm2; */
679c4762a1bSJed Brown /* IDi = (-QD*Vr + PD*Vi)/Vm2; */
680c4762a1bSJed Brown
681c4762a1bSJed Brown dIDr_dVr = (dPD_dVr * Vr + dQD_dVr * Vi + PD) / Vm2 - ((PD * Vr + QD * Vi) * 2 * Vr) / Vm4;
682c4762a1bSJed Brown dIDr_dVi = (dPD_dVi * Vr + dQD_dVi * Vi + QD) / Vm2 - ((PD * Vr + QD * Vi) * 2 * Vi) / Vm4;
683c4762a1bSJed Brown
684c4762a1bSJed Brown dIDi_dVr = (-dQD_dVr * Vr + dPD_dVr * Vi - QD) / Vm2 - ((-QD * Vr + PD * Vi) * 2 * Vr) / Vm4;
685c4762a1bSJed Brown dIDi_dVi = (-dQD_dVi * Vr + dPD_dVi * Vi + PD) / Vm2 - ((-QD * Vr + PD * Vi) * 2 * Vi) / Vm4;
686c4762a1bSJed Brown
687c4762a1bSJed Brown /* fnet[2*lbus[i]] += IDi; */
688c4762a1bSJed Brown row[0] = net_start + 2 * lbus[i];
6899371c9d4SSatish Balay col[0] = net_start + 2 * lbus[i];
6909371c9d4SSatish Balay col[1] = net_start + 2 * lbus[i] + 1;
6919371c9d4SSatish Balay val[0] = dIDi_dVr;
6929371c9d4SSatish Balay val[1] = dIDi_dVi;
6939566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, row, 2, col, val, ADD_VALUES));
694c4762a1bSJed Brown /* fnet[2*lbus[i]+1] += IDr; */
695c4762a1bSJed Brown row[0] = net_start + 2 * lbus[i] + 1;
6969371c9d4SSatish Balay col[0] = net_start + 2 * lbus[i];
6979371c9d4SSatish Balay col[1] = net_start + 2 * lbus[i] + 1;
6989371c9d4SSatish Balay val[0] = dIDr_dVr;
6999371c9d4SSatish Balay val[1] = dIDr_dVi;
7009566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, row, 2, col, val, ADD_VALUES));
701c4762a1bSJed Brown }
7029566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(user->V0, &v0));
703c4762a1bSJed Brown
7049566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Xgen, &xgen));
7059566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Xnet, &xnet));
706c4762a1bSJed Brown
7079566063dSJacob Faibussowitsch PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid, &Xgen, &Xnet));
708c4762a1bSJed Brown
7099566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
7109566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
7113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
712c4762a1bSJed Brown }
713c4762a1bSJed Brown
714c4762a1bSJed Brown /*
715c4762a1bSJed Brown J = [I, 0
716c4762a1bSJed Brown dg_dx, dg_dy]
717c4762a1bSJed Brown */
AlgJacobian(SNES snes,Vec X,Mat A,Mat B,PetscCtx ctx)718*2a8381b2SBarry Smith PetscErrorCode AlgJacobian(SNES snes, Vec X, Mat A, Mat B, PetscCtx ctx)
719d71ae5a4SJacob Faibussowitsch {
720c4762a1bSJed Brown Userctx *user = (Userctx *)ctx;
721c4762a1bSJed Brown
722c4762a1bSJed Brown PetscFunctionBegin;
7239566063dSJacob Faibussowitsch PetscCall(ResidualJacobian(snes, X, A, B, ctx));
7249566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_KEEP_NONZERO_PATTERN, PETSC_TRUE));
7259566063dSJacob Faibussowitsch PetscCall(MatZeroRowsIS(A, user->is_diff, 1.0, NULL, NULL));
7263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
727c4762a1bSJed Brown }
728c4762a1bSJed Brown
729c4762a1bSJed Brown /*
730c4762a1bSJed Brown J = [a*I-df_dx, -df_dy
731c4762a1bSJed Brown dg_dx, dg_dy]
732c4762a1bSJed Brown */
733c4762a1bSJed Brown
IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,Userctx * user)734d71ae5a4SJacob Faibussowitsch PetscErrorCode IJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat A, Mat B, Userctx *user)
735d71ae5a4SJacob Faibussowitsch {
736c4762a1bSJed Brown SNES snes;
737c4762a1bSJed Brown PetscScalar atmp = (PetscScalar)a;
738c4762a1bSJed Brown PetscInt i, row;
739c4762a1bSJed Brown
740c4762a1bSJed Brown PetscFunctionBegin;
741c4762a1bSJed Brown user->t = t;
742c4762a1bSJed Brown
7439566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes));
7449566063dSJacob Faibussowitsch PetscCall(ResidualJacobian(snes, X, A, B, user));
745c4762a1bSJed Brown for (i = 0; i < ngen; i++) {
746c4762a1bSJed Brown row = 9 * i;
7479566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES));
748c4762a1bSJed Brown row = 9 * i + 1;
7499566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES));
750c4762a1bSJed Brown row = 9 * i + 2;
7519566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES));
752c4762a1bSJed Brown row = 9 * i + 3;
7539566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES));
754c4762a1bSJed Brown row = 9 * i + 6;
7559566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES));
756c4762a1bSJed Brown row = 9 * i + 7;
7579566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES));
758c4762a1bSJed Brown row = 9 * i + 8;
7599566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES));
760c4762a1bSJed Brown }
7619566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
7629566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
7633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
764c4762a1bSJed Brown }
765c4762a1bSJed Brown
main(int argc,char ** argv)766d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
767d71ae5a4SJacob Faibussowitsch {
768c4762a1bSJed Brown TS ts;
769c4762a1bSJed Brown SNES snes_alg;
770c4762a1bSJed Brown PetscMPIInt size;
771c4762a1bSJed Brown Userctx user;
772c4762a1bSJed Brown PetscViewer Xview, Ybusview;
773c4762a1bSJed Brown Vec X;
774c4762a1bSJed Brown Mat J;
775c4762a1bSJed Brown PetscInt i;
776c4762a1bSJed Brown /* sensitivity context */
777c4762a1bSJed Brown PetscScalar *y_ptr;
778c4762a1bSJed Brown Vec lambda[1];
779c4762a1bSJed Brown PetscInt *idx2;
780c4762a1bSJed Brown Vec Xdot;
781c4762a1bSJed Brown Vec F_alg;
782c4762a1bSJed Brown PetscInt row_loc, col_loc;
783c4762a1bSJed Brown PetscScalar val;
784c4762a1bSJed Brown
785327415f7SBarry Smith PetscFunctionBeginUser;
7869566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, "petscoptions", help));
7879566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
7883c633725SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs");
789c4762a1bSJed Brown
790c4762a1bSJed Brown user.neqs_gen = 9 * ngen; /* # eqs. for generator subsystem */
791c4762a1bSJed Brown user.neqs_net = 2 * nbus; /* # eqs. for network subsystem */
792c4762a1bSJed Brown user.neqs_pgrid = user.neqs_gen + user.neqs_net;
793c4762a1bSJed Brown
794c4762a1bSJed Brown /* Create indices for differential and algebraic equations */
7959566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(7 * ngen, &idx2));
796c4762a1bSJed Brown for (i = 0; i < ngen; i++) {
7979371c9d4SSatish Balay idx2[7 * i] = 9 * i;
7989371c9d4SSatish Balay idx2[7 * i + 1] = 9 * i + 1;
7999371c9d4SSatish Balay idx2[7 * i + 2] = 9 * i + 2;
8009371c9d4SSatish Balay idx2[7 * i + 3] = 9 * i + 3;
8019371c9d4SSatish Balay idx2[7 * i + 4] = 9 * i + 6;
8029371c9d4SSatish Balay idx2[7 * i + 5] = 9 * i + 7;
8039371c9d4SSatish Balay idx2[7 * i + 6] = 9 * i + 8;
804c4762a1bSJed Brown }
8059566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, 7 * ngen, idx2, PETSC_COPY_VALUES, &user.is_diff));
8069566063dSJacob Faibussowitsch PetscCall(ISComplement(user.is_diff, 0, user.neqs_pgrid, &user.is_alg));
8079566063dSJacob Faibussowitsch PetscCall(PetscFree(idx2));
808c4762a1bSJed Brown
809c4762a1bSJed Brown /* Read initial voltage vector and Ybus */
8109566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "X.bin", FILE_MODE_READ, &Xview));
8119566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "Ybus.bin", FILE_MODE_READ, &Ybusview));
812c4762a1bSJed Brown
8139566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &user.V0));
8149566063dSJacob Faibussowitsch PetscCall(VecSetSizes(user.V0, PETSC_DECIDE, user.neqs_net));
8159566063dSJacob Faibussowitsch PetscCall(VecLoad(user.V0, Xview));
816c4762a1bSJed Brown
8179566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &user.Ybus));
8189566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user.Ybus, PETSC_DECIDE, PETSC_DECIDE, user.neqs_net, user.neqs_net));
8199566063dSJacob Faibussowitsch PetscCall(MatSetType(user.Ybus, MATBAIJ));
8209566063dSJacob Faibussowitsch /* PetscCall(MatSetBlockSize(user.Ybus,2)); */
8219566063dSJacob Faibussowitsch PetscCall(MatLoad(user.Ybus, Ybusview));
822c4762a1bSJed Brown
823c4762a1bSJed Brown /* Set run time options */
824d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Transient stability fault options", "");
825c4762a1bSJed Brown {
826c4762a1bSJed Brown user.tfaulton = 1.0;
827c4762a1bSJed Brown user.tfaultoff = 1.2;
828c4762a1bSJed Brown user.Rfault = 0.0001;
829c4762a1bSJed Brown user.faultbus = 8;
8309566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tfaulton", "", "", user.tfaulton, &user.tfaulton, NULL));
8319566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tfaultoff", "", "", user.tfaultoff, &user.tfaultoff, NULL));
8329566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-faultbus", "", "", user.faultbus, &user.faultbus, NULL));
833c4762a1bSJed Brown user.t0 = 0.0;
834c4762a1bSJed Brown user.tmax = 5.0;
8359566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-t0", "", "", user.t0, &user.t0, NULL));
8369566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tmax", "", "", user.tmax, &user.tmax, NULL));
837c4762a1bSJed Brown }
838d0609cedSBarry Smith PetscOptionsEnd();
839c4762a1bSJed Brown
8409566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&Xview));
8419566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&Ybusview));
842c4762a1bSJed Brown
843c4762a1bSJed Brown /* Create DMs for generator and network subsystems */
8449566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, user.neqs_gen, 1, 1, NULL, &user.dmgen));
8459566063dSJacob Faibussowitsch PetscCall(DMSetOptionsPrefix(user.dmgen, "dmgen_"));
8469566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(user.dmgen));
8479566063dSJacob Faibussowitsch PetscCall(DMSetUp(user.dmgen));
8489566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, user.neqs_net, 1, 1, NULL, &user.dmnet));
8499566063dSJacob Faibussowitsch PetscCall(DMSetOptionsPrefix(user.dmnet, "dmnet_"));
8509566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(user.dmnet));
8519566063dSJacob Faibussowitsch PetscCall(DMSetUp(user.dmnet));
852c4762a1bSJed Brown /* Create a composite DM packer and add the two DMs */
8539566063dSJacob Faibussowitsch PetscCall(DMCompositeCreate(PETSC_COMM_WORLD, &user.dmpgrid));
8549566063dSJacob Faibussowitsch PetscCall(DMSetOptionsPrefix(user.dmpgrid, "pgrid_"));
8559566063dSJacob Faibussowitsch PetscCall(DMCompositeAddDM(user.dmpgrid, user.dmgen));
8569566063dSJacob Faibussowitsch PetscCall(DMCompositeAddDM(user.dmpgrid, user.dmnet));
857c4762a1bSJed Brown
8589566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(user.dmpgrid, &X));
859c4762a1bSJed Brown
8609566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &J));
8619566063dSJacob Faibussowitsch PetscCall(MatSetSizes(J, PETSC_DECIDE, PETSC_DECIDE, user.neqs_pgrid, user.neqs_pgrid));
8629566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(J));
8639566063dSJacob Faibussowitsch PetscCall(PreallocateJacobian(J, &user));
864c4762a1bSJed Brown
865c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
866c4762a1bSJed Brown Create timestepping solver context
867c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
8689566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
8699566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
8709566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSCN));
8718434afd1SBarry Smith PetscCall(TSSetIFunction(ts, NULL, (TSIFunctionFn *)IFunction, &user));
8728434afd1SBarry Smith PetscCall(TSSetIJacobian(ts, J, J, (TSIJacobianFn *)IJacobian, &user));
8739566063dSJacob Faibussowitsch PetscCall(TSSetApplicationContext(ts, &user));
874c4762a1bSJed Brown
875c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
876c4762a1bSJed Brown Set initial conditions
877c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
8789566063dSJacob Faibussowitsch PetscCall(SetInitialGuess(X, &user));
879c4762a1bSJed Brown /* Just to set up the Jacobian structure */
8809566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X, &Xdot));
8819566063dSJacob Faibussowitsch PetscCall(IJacobian(ts, 0.0, X, Xdot, 0.0, J, J, &user));
8829566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Xdot));
883c4762a1bSJed Brown
884c4762a1bSJed Brown /*
885c4762a1bSJed Brown Save trajectory of solution so that TSAdjointSolve() may be used
886c4762a1bSJed Brown */
8879566063dSJacob Faibussowitsch PetscCall(TSSetSaveTrajectory(ts));
888c4762a1bSJed Brown
8899566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, user.tmax));
8909566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, 0.01));
8919566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
8929566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts));
893c4762a1bSJed Brown
894c4762a1bSJed Brown user.alg_flg = PETSC_FALSE;
895c4762a1bSJed Brown /* Prefault period */
8969566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, X));
897c4762a1bSJed Brown
898c4762a1bSJed Brown /* Create the nonlinear solver for solving the algebraic system */
899c4762a1bSJed Brown /* Note that although the algebraic system needs to be solved only for
900c4762a1bSJed Brown Idq and V, we reuse the entire system including xgen. The xgen
901c4762a1bSJed Brown variables are held constant by setting their residuals to 0 and
902c4762a1bSJed Brown putting a 1 on the Jacobian diagonal for xgen rows
903c4762a1bSJed Brown */
9049566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X, &F_alg));
9059566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes_alg));
9069566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes_alg, F_alg, AlgFunction, &user));
9079566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(J));
9089566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes_alg, J, J, AlgJacobian, &user));
9099566063dSJacob Faibussowitsch PetscCall(SNESSetOptionsPrefix(snes_alg, "alg_"));
9109566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes_alg));
911c4762a1bSJed Brown
912c4762a1bSJed Brown /* Apply disturbance - resistive fault at user.faultbus */
913c4762a1bSJed Brown /* This is done by adding shunt conductance to the diagonal location
914c4762a1bSJed Brown in the Ybus matrix */
9159371c9d4SSatish Balay row_loc = 2 * user.faultbus;
9169371c9d4SSatish Balay col_loc = 2 * user.faultbus + 1; /* Location for G */
917c4762a1bSJed Brown val = 1 / user.Rfault;
9189566063dSJacob Faibussowitsch PetscCall(MatSetValues(user.Ybus, 1, &row_loc, 1, &col_loc, &val, ADD_VALUES));
9199371c9d4SSatish Balay row_loc = 2 * user.faultbus + 1;
9209371c9d4SSatish Balay col_loc = 2 * user.faultbus; /* Location for G */
921c4762a1bSJed Brown val = 1 / user.Rfault;
9229566063dSJacob Faibussowitsch PetscCall(MatSetValues(user.Ybus, 1, &row_loc, 1, &col_loc, &val, ADD_VALUES));
923c4762a1bSJed Brown
9249566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user.Ybus, MAT_FINAL_ASSEMBLY));
9259566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user.Ybus, MAT_FINAL_ASSEMBLY));
926c4762a1bSJed Brown
927c4762a1bSJed Brown user.alg_flg = PETSC_TRUE;
928c4762a1bSJed Brown /* Solve the algebraic equations */
9299566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes_alg, NULL, X));
930c4762a1bSJed Brown
931c4762a1bSJed Brown /* Disturbance period */
932c4762a1bSJed Brown user.alg_flg = PETSC_FALSE;
9339566063dSJacob Faibussowitsch PetscCall(TSSetTime(ts, user.tfaulton));
9349566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, user.tfaultoff));
9359566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, X));
936c4762a1bSJed Brown
937c4762a1bSJed Brown /* Remove the fault */
9389371c9d4SSatish Balay row_loc = 2 * user.faultbus;
9399371c9d4SSatish Balay col_loc = 2 * user.faultbus + 1;
940c4762a1bSJed Brown val = -1 / user.Rfault;
9419566063dSJacob Faibussowitsch PetscCall(MatSetValues(user.Ybus, 1, &row_loc, 1, &col_loc, &val, ADD_VALUES));
9429371c9d4SSatish Balay row_loc = 2 * user.faultbus + 1;
9439371c9d4SSatish Balay col_loc = 2 * user.faultbus;
944c4762a1bSJed Brown val = -1 / user.Rfault;
9459566063dSJacob Faibussowitsch PetscCall(MatSetValues(user.Ybus, 1, &row_loc, 1, &col_loc, &val, ADD_VALUES));
946c4762a1bSJed Brown
9479566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(user.Ybus, MAT_FINAL_ASSEMBLY));
9489566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(user.Ybus, MAT_FINAL_ASSEMBLY));
949c4762a1bSJed Brown
9509566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(J));
951c4762a1bSJed Brown
952c4762a1bSJed Brown user.alg_flg = PETSC_TRUE;
953c4762a1bSJed Brown
954c4762a1bSJed Brown /* Solve the algebraic equations */
9559566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes_alg, NULL, X));
956c4762a1bSJed Brown
957c4762a1bSJed Brown /* Post-disturbance period */
958c4762a1bSJed Brown user.alg_flg = PETSC_TRUE;
9599566063dSJacob Faibussowitsch PetscCall(TSSetTime(ts, user.tfaultoff));
9609566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, user.tmax));
9619566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, X));
962c4762a1bSJed Brown
963c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
964c4762a1bSJed Brown Adjoint model starts here
965c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
9669566063dSJacob Faibussowitsch PetscCall(TSSetPostStep(ts, NULL));
9679566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(J, &lambda[0], NULL));
968c4762a1bSJed Brown /* Set initial conditions for the adjoint integration */
9699566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(lambda[0]));
9709566063dSJacob Faibussowitsch PetscCall(VecGetArray(lambda[0], &y_ptr));
971c4762a1bSJed Brown y_ptr[0] = 1.0;
9729566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(lambda[0], &y_ptr));
9739566063dSJacob Faibussowitsch PetscCall(TSSetCostGradients(ts, 1, lambda, NULL));
974c4762a1bSJed Brown
9759566063dSJacob Faibussowitsch PetscCall(TSAdjointSolve(ts));
976c4762a1bSJed Brown
9779566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n sensitivity wrt initial conditions: \n"));
9789566063dSJacob Faibussowitsch PetscCall(VecView(lambda[0], PETSC_VIEWER_STDOUT_WORLD));
9799566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lambda[0]));
980c4762a1bSJed Brown
9819566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes_alg));
9829566063dSJacob Faibussowitsch PetscCall(VecDestroy(&F_alg));
9839566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J));
9849566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user.Ybus));
9859566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X));
9869566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user.V0));
9879566063dSJacob Faibussowitsch PetscCall(DMDestroy(&user.dmgen));
9889566063dSJacob Faibussowitsch PetscCall(DMDestroy(&user.dmnet));
9899566063dSJacob Faibussowitsch PetscCall(DMDestroy(&user.dmpgrid));
9909566063dSJacob Faibussowitsch PetscCall(ISDestroy(&user.is_diff));
9919566063dSJacob Faibussowitsch PetscCall(ISDestroy(&user.is_alg));
9929566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts));
9939566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
994b122ec5aSJacob Faibussowitsch return 0;
995c4762a1bSJed Brown }
996c4762a1bSJed Brown
997c4762a1bSJed Brown /*TEST
998c4762a1bSJed Brown
999c4762a1bSJed Brown build:
1000dfd57a17SPierre Jolivet requires: double !complex !defined(PETSC_USE_64BIT_INDICES)
1001c4762a1bSJed Brown
1002c4762a1bSJed Brown test:
1003c4762a1bSJed Brown args: -viewer_binary_skip_info
1004c4762a1bSJed Brown localrunfiles: petscoptions X.bin Ybus.bin
1005c4762a1bSJed Brown
100638bc3bfcSIvan Yashchuk test:
100738bc3bfcSIvan Yashchuk suffix: 2
100838bc3bfcSIvan Yashchuk args: -viewer_binary_skip_info -ts_type beuler
100938bc3bfcSIvan Yashchuk localrunfiles: petscoptions X.bin Ybus.bin
101038bc3bfcSIvan Yashchuk
1011c4762a1bSJed Brown TEST*/
1012