1c4762a1bSJed Brown static char help[] = "Using finite difference for the problem in ex9busopt.c \n\n";
2c4762a1bSJed Brown
3c4762a1bSJed Brown /*
4c4762a1bSJed Brown Use finite difference approximations to solve the same optimization problem as in ex9busopt.c.
5c4762a1bSJed Brown */
6c4762a1bSJed Brown
7c4762a1bSJed Brown #include <petsctao.h>
8c4762a1bSJed Brown #include <petscts.h>
9c4762a1bSJed Brown #include <petscdm.h>
10c4762a1bSJed Brown #include <petscdmda.h>
11c4762a1bSJed Brown #include <petscdmcomposite.h>
12c4762a1bSJed Brown
13c4762a1bSJed Brown PetscErrorCode FormFunction(Tao, Vec, PetscReal *, void *);
14c4762a1bSJed Brown
15c4762a1bSJed Brown #define freq 60
16c4762a1bSJed Brown #define w_s (2 * PETSC_PI * freq)
17c4762a1bSJed Brown
18c4762a1bSJed Brown /* Sizes and indices */
19c4762a1bSJed Brown const PetscInt nbus = 9; /* Number of network buses */
20c4762a1bSJed Brown const PetscInt ngen = 3; /* Number of generators */
21c4762a1bSJed Brown const PetscInt nload = 3; /* Number of loads */
22c4762a1bSJed Brown const PetscInt gbus[3] = {0, 1, 2}; /* Buses at which generators are incident */
23c4762a1bSJed Brown const PetscInt lbus[3] = {4, 5, 7}; /* Buses at which loads are incident */
24c4762a1bSJed Brown
25c4762a1bSJed Brown /* Generator real and reactive powers (found via loadflow) */
26c4762a1bSJed Brown PetscScalar PG[3] = {0.69, 1.59, 0.69};
27c4762a1bSJed Brown /* PetscScalar PG[3] = {0.716786142395021,1.630000000000000,0.850000000000000};*/
28c4762a1bSJed Brown const PetscScalar QG[3] = {0.270702180178785, 0.066120127797275, -0.108402221791588};
29c4762a1bSJed Brown /* Generator constants */
30c4762a1bSJed Brown const PetscScalar H[3] = {23.64, 6.4, 3.01}; /* Inertia constant */
31c4762a1bSJed Brown const PetscScalar Rs[3] = {0.0, 0.0, 0.0}; /* Stator Resistance */
32c4762a1bSJed Brown const PetscScalar Xd[3] = {0.146, 0.8958, 1.3125}; /* d-axis reactance */
33c4762a1bSJed Brown const PetscScalar Xdp[3] = {0.0608, 0.1198, 0.1813}; /* d-axis transient reactance */
34c4762a1bSJed 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 */
35c4762a1bSJed Brown const PetscScalar Xqp[3] = {0.0969, 0.1969, 0.25}; /* q-axis transient reactance */
36c4762a1bSJed Brown const PetscScalar Td0p[3] = {8.96, 6.0, 5.89}; /* d-axis open circuit time constant */
37c4762a1bSJed Brown const PetscScalar Tq0p[3] = {0.31, 0.535, 0.6}; /* q-axis open circuit time constant */
38c4762a1bSJed Brown PetscScalar M[3]; /* M = 2*H/w_s */
39c4762a1bSJed Brown PetscScalar D[3]; /* D = 0.1*M */
40c4762a1bSJed Brown
41c4762a1bSJed Brown PetscScalar TM[3]; /* Mechanical Torque */
42c4762a1bSJed Brown /* Exciter system constants */
43c4762a1bSJed Brown const PetscScalar KA[3] = {20.0, 20.0, 20.0}; /* Voltage regulartor gain constant */
44c4762a1bSJed Brown const PetscScalar TA[3] = {0.2, 0.2, 0.2}; /* Voltage regulator time constant */
45c4762a1bSJed Brown const PetscScalar KE[3] = {1.0, 1.0, 1.0}; /* Exciter gain constant */
46c4762a1bSJed Brown const PetscScalar TE[3] = {0.314, 0.314, 0.314}; /* Exciter time constant */
47c4762a1bSJed Brown const PetscScalar KF[3] = {0.063, 0.063, 0.063}; /* Feedback stabilizer gain constant */
48c4762a1bSJed Brown const PetscScalar TF[3] = {0.35, 0.35, 0.35}; /* Feedback stabilizer time constant */
49c4762a1bSJed Brown const PetscScalar k1[3] = {0.0039, 0.0039, 0.0039};
50c4762a1bSJed Brown const PetscScalar k2[3] = {1.555, 1.555, 1.555}; /* k1 and k2 for calculating the saturation function SE = k1*exp(k2*Efd) */
51c4762a1bSJed Brown
52c4762a1bSJed Brown PetscScalar Vref[3];
53c4762a1bSJed Brown /* Load constants
54c4762a1bSJed Brown We use a composite load model that describes the load and reactive powers at each time instant as follows
55c4762a1bSJed Brown P(t) = \sum\limits_{i=0}^ld_nsegsp \ld_alphap_i*P_D0(\frac{V_m(t)}{V_m0})^\ld_betap_i
56c4762a1bSJed Brown Q(t) = \sum\limits_{i=0}^ld_nsegsq \ld_alphaq_i*Q_D0(\frac{V_m(t)}{V_m0})^\ld_betaq_i
57c4762a1bSJed Brown where
58c4762a1bSJed Brown ld_nsegsp,ld_nsegsq - Number of individual load models for real and reactive power loads
59c4762a1bSJed Brown ld_alphap,ld_alphap - Percentage contribution (weights) or loads
60c4762a1bSJed Brown P_D0 - Real power load
61c4762a1bSJed Brown Q_D0 - Reactive power load
62c4762a1bSJed Brown V_m(t) - Voltage magnitude at time t
63c4762a1bSJed Brown V_m0 - Voltage magnitude at t = 0
64c4762a1bSJed Brown ld_betap, ld_betaq - exponents describing the load model for real and reactive part
65c4762a1bSJed Brown
66c4762a1bSJed Brown Note: All loads have the same characteristic currently.
67c4762a1bSJed Brown */
68c4762a1bSJed Brown const PetscScalar PD0[3] = {1.25, 0.9, 1.0};
69c4762a1bSJed Brown const PetscScalar QD0[3] = {0.5, 0.3, 0.35};
70c4762a1bSJed Brown const PetscInt ld_nsegsp[3] = {3, 3, 3};
71c4762a1bSJed Brown const PetscScalar ld_alphap[3] = {1.0, 0.0, 0.0};
72c4762a1bSJed Brown const PetscScalar ld_betap[3] = {2.0, 1.0, 0.0};
73c4762a1bSJed Brown const PetscInt ld_nsegsq[3] = {3, 3, 3};
74c4762a1bSJed Brown const PetscScalar ld_alphaq[3] = {1.0, 0.0, 0.0};
75c4762a1bSJed Brown const PetscScalar ld_betaq[3] = {2.0, 1.0, 0.0};
76c4762a1bSJed Brown
77c4762a1bSJed Brown typedef struct {
78c4762a1bSJed Brown DM dmgen, dmnet; /* DMs to manage generator and network subsystem */
79c4762a1bSJed Brown DM dmpgrid; /* Composite DM to manage the entire power grid */
80c4762a1bSJed Brown Mat Ybus; /* Network admittance matrix */
81c4762a1bSJed Brown Vec V0; /* Initial voltage vector (Power flow solution) */
82c4762a1bSJed Brown PetscReal tfaulton, tfaultoff; /* Fault on and off times */
83c4762a1bSJed Brown PetscInt faultbus; /* Fault bus */
84c4762a1bSJed Brown PetscScalar Rfault;
85c4762a1bSJed Brown PetscReal t0, tmax;
86c4762a1bSJed Brown PetscInt neqs_gen, neqs_net, neqs_pgrid;
87c4762a1bSJed Brown Mat Sol; /* Matrix to save solution at each time step */
88c4762a1bSJed Brown PetscInt stepnum;
89c4762a1bSJed Brown PetscBool alg_flg;
90c4762a1bSJed Brown PetscReal t;
91c4762a1bSJed Brown IS is_diff; /* indices for differential equations */
92c4762a1bSJed Brown IS is_alg; /* indices for algebraic equations */
93c4762a1bSJed Brown PetscReal freq_u, freq_l; /* upper and lower frequency limit */
94c4762a1bSJed Brown PetscInt pow; /* power coefficient used in the cost function */
95c4762a1bSJed Brown Vec vec_q;
96c4762a1bSJed Brown } Userctx;
97c4762a1bSJed Brown
98c4762a1bSJed 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)99d71ae5a4SJacob Faibussowitsch PetscErrorCode dq2ri(PetscScalar Fd, PetscScalar Fq, PetscScalar delta, PetscScalar *Fr, PetscScalar *Fi)
100d71ae5a4SJacob Faibussowitsch {
101c4762a1bSJed Brown PetscFunctionBegin;
102c4762a1bSJed Brown *Fr = Fd * PetscSinScalar(delta) + Fq * PetscCosScalar(delta);
103c4762a1bSJed Brown *Fi = -Fd * PetscCosScalar(delta) + Fq * PetscSinScalar(delta);
1043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
105c4762a1bSJed Brown }
106c4762a1bSJed Brown
107c4762a1bSJed 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)108d71ae5a4SJacob Faibussowitsch PetscErrorCode ri2dq(PetscScalar Fr, PetscScalar Fi, PetscScalar delta, PetscScalar *Fd, PetscScalar *Fq)
109d71ae5a4SJacob Faibussowitsch {
110c4762a1bSJed Brown PetscFunctionBegin;
111c4762a1bSJed Brown *Fd = Fr * PetscSinScalar(delta) - Fi * PetscCosScalar(delta);
112c4762a1bSJed Brown *Fq = Fr * PetscCosScalar(delta) + Fi * PetscSinScalar(delta);
1133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
114c4762a1bSJed Brown }
115c4762a1bSJed Brown
SetInitialGuess(Vec X,Userctx * user)116d71ae5a4SJacob Faibussowitsch PetscErrorCode SetInitialGuess(Vec X, Userctx *user)
117d71ae5a4SJacob Faibussowitsch {
118c4762a1bSJed Brown Vec Xgen, Xnet;
119c4762a1bSJed Brown PetscScalar *xgen, *xnet;
120c4762a1bSJed Brown PetscInt i, idx = 0;
121c4762a1bSJed Brown PetscScalar Vr, Vi, IGr, IGi, Vm, Vm2;
122c4762a1bSJed Brown PetscScalar Eqp, Edp, delta;
123c4762a1bSJed Brown PetscScalar Efd, RF, VR; /* Exciter variables */
124c4762a1bSJed Brown PetscScalar Id, Iq; /* Generator dq axis currents */
125c4762a1bSJed Brown PetscScalar theta, Vd, Vq, SE;
126c4762a1bSJed Brown
127c4762a1bSJed Brown PetscFunctionBegin;
1289371c9d4SSatish Balay M[0] = 2 * H[0] / w_s;
1299371c9d4SSatish Balay M[1] = 2 * H[1] / w_s;
1309371c9d4SSatish Balay M[2] = 2 * H[2] / w_s;
1319371c9d4SSatish Balay D[0] = 0.1 * M[0];
1329371c9d4SSatish Balay D[1] = 0.1 * M[1];
1339371c9d4SSatish Balay D[2] = 0.1 * M[2];
134c4762a1bSJed Brown
1359566063dSJacob Faibussowitsch PetscCall(DMCompositeGetLocalVectors(user->dmpgrid, &Xgen, &Xnet));
136c4762a1bSJed Brown
137c4762a1bSJed Brown /* Network subsystem initialization */
1389566063dSJacob Faibussowitsch PetscCall(VecCopy(user->V0, Xnet));
139c4762a1bSJed Brown
140c4762a1bSJed Brown /* Generator subsystem initialization */
1419566063dSJacob Faibussowitsch PetscCall(VecGetArray(Xgen, &xgen));
1429566063dSJacob Faibussowitsch PetscCall(VecGetArray(Xnet, &xnet));
143c4762a1bSJed Brown
144c4762a1bSJed Brown for (i = 0; i < ngen; i++) {
145c4762a1bSJed Brown Vr = xnet[2 * gbus[i]]; /* Real part of generator terminal voltage */
146c4762a1bSJed Brown Vi = xnet[2 * gbus[i] + 1]; /* Imaginary part of the generator terminal voltage */
1479371c9d4SSatish Balay Vm = PetscSqrtScalar(Vr * Vr + Vi * Vi);
1489371c9d4SSatish Balay Vm2 = Vm * Vm;
149c4762a1bSJed Brown IGr = (Vr * PG[i] + Vi * QG[i]) / Vm2;
150c4762a1bSJed Brown IGi = (Vi * PG[i] - Vr * QG[i]) / Vm2;
151c4762a1bSJed Brown
152c4762a1bSJed Brown delta = PetscAtan2Real(Vi + Xq[i] * IGr, Vr - Xq[i] * IGi); /* Machine angle */
153c4762a1bSJed Brown
154c4762a1bSJed Brown theta = PETSC_PI / 2.0 - delta;
155c4762a1bSJed Brown
156c4762a1bSJed Brown Id = IGr * PetscCosScalar(theta) - IGi * PetscSinScalar(theta); /* d-axis stator current */
157c4762a1bSJed Brown Iq = IGr * PetscSinScalar(theta) + IGi * PetscCosScalar(theta); /* q-axis stator current */
158c4762a1bSJed Brown
159c4762a1bSJed Brown Vd = Vr * PetscCosScalar(theta) - Vi * PetscSinScalar(theta);
160c4762a1bSJed Brown Vq = Vr * PetscSinScalar(theta) + Vi * PetscCosScalar(theta);
161c4762a1bSJed Brown
162c4762a1bSJed Brown Edp = Vd + Rs[i] * Id - Xqp[i] * Iq; /* d-axis transient EMF */
163c4762a1bSJed Brown Eqp = Vq + Rs[i] * Iq + Xdp[i] * Id; /* q-axis transient EMF */
164c4762a1bSJed Brown
165c4762a1bSJed Brown TM[i] = PG[i];
166c4762a1bSJed Brown
167c4762a1bSJed Brown /* The generator variables are ordered as [Eqp,Edp,delta,w,Id,Iq] */
168c4762a1bSJed Brown xgen[idx] = Eqp;
169c4762a1bSJed Brown xgen[idx + 1] = Edp;
170c4762a1bSJed Brown xgen[idx + 2] = delta;
171c4762a1bSJed Brown xgen[idx + 3] = w_s;
172c4762a1bSJed Brown
173c4762a1bSJed Brown idx = idx + 4;
174c4762a1bSJed Brown
175c4762a1bSJed Brown xgen[idx] = Id;
176c4762a1bSJed Brown xgen[idx + 1] = Iq;
177c4762a1bSJed Brown
178c4762a1bSJed Brown idx = idx + 2;
179c4762a1bSJed Brown
180c4762a1bSJed Brown /* Exciter */
181c4762a1bSJed Brown Efd = Eqp + (Xd[i] - Xdp[i]) * Id;
182c4762a1bSJed Brown SE = k1[i] * PetscExpScalar(k2[i] * Efd);
183c4762a1bSJed Brown VR = KE[i] * Efd + SE;
184c4762a1bSJed Brown RF = KF[i] * Efd / TF[i];
185c4762a1bSJed Brown
186c4762a1bSJed Brown xgen[idx] = Efd;
187c4762a1bSJed Brown xgen[idx + 1] = RF;
188c4762a1bSJed Brown xgen[idx + 2] = VR;
189c4762a1bSJed Brown
190c4762a1bSJed Brown Vref[i] = Vm + (VR / KA[i]);
191c4762a1bSJed Brown
192c4762a1bSJed Brown idx = idx + 3;
193c4762a1bSJed Brown }
194c4762a1bSJed Brown
1959566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Xgen, &xgen));
1969566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Xnet, &xnet));
197c4762a1bSJed Brown
1989566063dSJacob Faibussowitsch /* PetscCall(VecView(Xgen,0)); */
1999566063dSJacob Faibussowitsch PetscCall(DMCompositeGather(user->dmpgrid, INSERT_VALUES, X, Xgen, Xnet));
2009566063dSJacob Faibussowitsch PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid, &Xgen, &Xnet));
2013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
202c4762a1bSJed Brown }
203c4762a1bSJed Brown
204c4762a1bSJed Brown /* Computes F = [-f(x,y);g(x,y)] */
ResidualFunction(SNES snes,Vec X,Vec F,Userctx * user)205d71ae5a4SJacob Faibussowitsch PetscErrorCode ResidualFunction(SNES snes, Vec X, Vec F, Userctx *user)
206d71ae5a4SJacob Faibussowitsch {
207c4762a1bSJed Brown Vec Xgen, Xnet, Fgen, Fnet;
208c4762a1bSJed Brown PetscScalar *xgen, *xnet, *fgen, *fnet;
209c4762a1bSJed Brown PetscInt i, idx = 0;
210c4762a1bSJed Brown PetscScalar Vr, Vi, Vm, Vm2;
211c4762a1bSJed Brown PetscScalar Eqp, Edp, delta, w; /* Generator variables */
212c4762a1bSJed Brown PetscScalar Efd, RF, VR; /* Exciter variables */
213c4762a1bSJed Brown PetscScalar Id, Iq; /* Generator dq axis currents */
214c4762a1bSJed Brown PetscScalar Vd, Vq, SE;
215c4762a1bSJed Brown PetscScalar IGr, IGi, IDr, IDi;
216c4762a1bSJed Brown PetscScalar Zdq_inv[4], det;
217c4762a1bSJed Brown PetscScalar PD, QD, Vm0, *v0;
218c4762a1bSJed Brown PetscInt k;
219c4762a1bSJed Brown
220c4762a1bSJed Brown PetscFunctionBegin;
2219566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F));
2229566063dSJacob Faibussowitsch PetscCall(DMCompositeGetLocalVectors(user->dmpgrid, &Xgen, &Xnet));
2239566063dSJacob Faibussowitsch PetscCall(DMCompositeGetLocalVectors(user->dmpgrid, &Fgen, &Fnet));
2249566063dSJacob Faibussowitsch PetscCall(DMCompositeScatter(user->dmpgrid, X, Xgen, Xnet));
2259566063dSJacob Faibussowitsch PetscCall(DMCompositeScatter(user->dmpgrid, F, Fgen, Fnet));
226c4762a1bSJed Brown
227c4762a1bSJed Brown /* Network current balance residual IG + Y*V + IL = 0. Only YV is added here.
228c4762a1bSJed Brown The generator current injection, IG, and load current injection, ID are added later
229c4762a1bSJed Brown */
230c4762a1bSJed Brown /* Note that the values in Ybus are stored assuming the imaginary current balance
231c4762a1bSJed Brown equation is ordered first followed by real current balance equation for each bus.
232c4762a1bSJed Brown Thus imaginary current contribution goes in location 2*i, and
233c4762a1bSJed Brown real current contribution in 2*i+1
234c4762a1bSJed Brown */
2359566063dSJacob Faibussowitsch PetscCall(MatMult(user->Ybus, Xnet, Fnet));
236c4762a1bSJed Brown
2379566063dSJacob Faibussowitsch PetscCall(VecGetArray(Xgen, &xgen));
2389566063dSJacob Faibussowitsch PetscCall(VecGetArray(Xnet, &xnet));
2399566063dSJacob Faibussowitsch PetscCall(VecGetArray(Fgen, &fgen));
2409566063dSJacob Faibussowitsch PetscCall(VecGetArray(Fnet, &fnet));
241c4762a1bSJed Brown
242c4762a1bSJed Brown /* Generator subsystem */
243c4762a1bSJed Brown for (i = 0; i < ngen; i++) {
244c4762a1bSJed Brown Eqp = xgen[idx];
245c4762a1bSJed Brown Edp = xgen[idx + 1];
246c4762a1bSJed Brown delta = xgen[idx + 2];
247c4762a1bSJed Brown w = xgen[idx + 3];
248c4762a1bSJed Brown Id = xgen[idx + 4];
249c4762a1bSJed Brown Iq = xgen[idx + 5];
250c4762a1bSJed Brown Efd = xgen[idx + 6];
251c4762a1bSJed Brown RF = xgen[idx + 7];
252c4762a1bSJed Brown VR = xgen[idx + 8];
253c4762a1bSJed Brown
254c4762a1bSJed Brown /* Generator differential equations */
255c4762a1bSJed Brown fgen[idx] = (Eqp + (Xd[i] - Xdp[i]) * Id - Efd) / Td0p[i];
256c4762a1bSJed Brown fgen[idx + 1] = (Edp - (Xq[i] - Xqp[i]) * Iq) / Tq0p[i];
257c4762a1bSJed Brown fgen[idx + 2] = -w + w_s;
258c4762a1bSJed Brown fgen[idx + 3] = (-TM[i] + Edp * Id + Eqp * Iq + (Xqp[i] - Xdp[i]) * Id * Iq + D[i] * (w - w_s)) / M[i];
259c4762a1bSJed Brown
260c4762a1bSJed Brown Vr = xnet[2 * gbus[i]]; /* Real part of generator terminal voltage */
261c4762a1bSJed Brown Vi = xnet[2 * gbus[i] + 1]; /* Imaginary part of the generator terminal voltage */
262c4762a1bSJed Brown
2639566063dSJacob Faibussowitsch PetscCall(ri2dq(Vr, Vi, delta, &Vd, &Vq));
264c4762a1bSJed Brown /* Algebraic equations for stator currents */
265c4762a1bSJed Brown det = Rs[i] * Rs[i] + Xdp[i] * Xqp[i];
266c4762a1bSJed Brown
267c4762a1bSJed Brown Zdq_inv[0] = Rs[i] / det;
268c4762a1bSJed Brown Zdq_inv[1] = Xqp[i] / det;
269c4762a1bSJed Brown Zdq_inv[2] = -Xdp[i] / det;
270c4762a1bSJed Brown Zdq_inv[3] = Rs[i] / det;
271c4762a1bSJed Brown
272c4762a1bSJed Brown fgen[idx + 4] = Zdq_inv[0] * (-Edp + Vd) + Zdq_inv[1] * (-Eqp + Vq) + Id;
273c4762a1bSJed Brown fgen[idx + 5] = Zdq_inv[2] * (-Edp + Vd) + Zdq_inv[3] * (-Eqp + Vq) + Iq;
274c4762a1bSJed Brown
275c4762a1bSJed Brown /* Add generator current injection to network */
2769566063dSJacob Faibussowitsch PetscCall(dq2ri(Id, Iq, delta, &IGr, &IGi));
277c4762a1bSJed Brown
278c4762a1bSJed Brown fnet[2 * gbus[i]] -= IGi;
279c4762a1bSJed Brown fnet[2 * gbus[i] + 1] -= IGr;
280c4762a1bSJed Brown
281c4762a1bSJed Brown Vm = PetscSqrtScalar(Vd * Vd + Vq * Vq);
282c4762a1bSJed Brown
283c4762a1bSJed Brown SE = k1[i] * PetscExpScalar(k2[i] * Efd);
284c4762a1bSJed Brown
285c4762a1bSJed Brown /* Exciter differential equations */
286c4762a1bSJed Brown fgen[idx + 6] = (KE[i] * Efd + SE - VR) / TE[i];
287c4762a1bSJed Brown fgen[idx + 7] = (RF - KF[i] * Efd / TF[i]) / TF[i];
288c4762a1bSJed Brown fgen[idx + 8] = (VR - KA[i] * RF + KA[i] * KF[i] * Efd / TF[i] - KA[i] * (Vref[i] - Vm)) / TA[i];
289c4762a1bSJed Brown
290c4762a1bSJed Brown idx = idx + 9;
291c4762a1bSJed Brown }
292c4762a1bSJed Brown
2939566063dSJacob Faibussowitsch PetscCall(VecGetArray(user->V0, &v0));
294c4762a1bSJed Brown for (i = 0; i < nload; i++) {
295c4762a1bSJed Brown Vr = xnet[2 * lbus[i]]; /* Real part of load bus voltage */
296c4762a1bSJed Brown Vi = xnet[2 * lbus[i] + 1]; /* Imaginary part of the load bus voltage */
2979371c9d4SSatish Balay Vm = PetscSqrtScalar(Vr * Vr + Vi * Vi);
2989371c9d4SSatish Balay Vm2 = Vm * Vm;
299c4762a1bSJed Brown Vm0 = PetscSqrtScalar(v0[2 * lbus[i]] * v0[2 * lbus[i]] + v0[2 * lbus[i] + 1] * v0[2 * lbus[i] + 1]);
300c4762a1bSJed Brown PD = QD = 0.0;
30157508eceSPierre Jolivet for (k = 0; k < ld_nsegsp[i]; k++) PD += ld_alphap[k] * PD0[i] * PetscPowScalar(Vm / Vm0, ld_betap[k]);
30257508eceSPierre Jolivet for (k = 0; k < ld_nsegsq[i]; k++) QD += ld_alphaq[k] * QD0[i] * PetscPowScalar(Vm / Vm0, ld_betaq[k]);
303c4762a1bSJed Brown
304c4762a1bSJed Brown /* Load currents */
305c4762a1bSJed Brown IDr = (PD * Vr + QD * Vi) / Vm2;
306c4762a1bSJed Brown IDi = (-QD * Vr + PD * Vi) / Vm2;
307c4762a1bSJed Brown
308c4762a1bSJed Brown fnet[2 * lbus[i]] += IDi;
309c4762a1bSJed Brown fnet[2 * lbus[i] + 1] += IDr;
310c4762a1bSJed Brown }
3119566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(user->V0, &v0));
312c4762a1bSJed Brown
3139566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Xgen, &xgen));
3149566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Xnet, &xnet));
3159566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Fgen, &fgen));
3169566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Fnet, &fnet));
317c4762a1bSJed Brown
3189566063dSJacob Faibussowitsch PetscCall(DMCompositeGather(user->dmpgrid, INSERT_VALUES, F, Fgen, Fnet));
3199566063dSJacob Faibussowitsch PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid, &Xgen, &Xnet));
3209566063dSJacob Faibussowitsch PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid, &Fgen, &Fnet));
3213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
322c4762a1bSJed Brown }
323c4762a1bSJed Brown
324c4762a1bSJed Brown /* \dot{x} - f(x,y)
325c4762a1bSJed Brown g(x,y) = 0
326c4762a1bSJed Brown */
IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,Userctx * user)327d71ae5a4SJacob Faibussowitsch PetscErrorCode IFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, Userctx *user)
328d71ae5a4SJacob Faibussowitsch {
329c4762a1bSJed Brown SNES snes;
330c4762a1bSJed Brown PetscScalar *f;
331c4762a1bSJed Brown const PetscScalar *xdot;
332c4762a1bSJed Brown PetscInt i;
333c4762a1bSJed Brown
334c4762a1bSJed Brown PetscFunctionBegin;
335c4762a1bSJed Brown user->t = t;
336c4762a1bSJed Brown
3379566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes));
3389566063dSJacob Faibussowitsch PetscCall(ResidualFunction(snes, X, F, user));
3399566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f));
3409566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Xdot, &xdot));
341c4762a1bSJed Brown for (i = 0; i < ngen; i++) {
342c4762a1bSJed Brown f[9 * i] += xdot[9 * i];
343c4762a1bSJed Brown f[9 * i + 1] += xdot[9 * i + 1];
344c4762a1bSJed Brown f[9 * i + 2] += xdot[9 * i + 2];
345c4762a1bSJed Brown f[9 * i + 3] += xdot[9 * i + 3];
346c4762a1bSJed Brown f[9 * i + 6] += xdot[9 * i + 6];
347c4762a1bSJed Brown f[9 * i + 7] += xdot[9 * i + 7];
348c4762a1bSJed Brown f[9 * i + 8] += xdot[9 * i + 8];
349c4762a1bSJed Brown }
3509566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f));
3519566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Xdot, &xdot));
3523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
353c4762a1bSJed Brown }
354c4762a1bSJed Brown
355c4762a1bSJed Brown /* This function is used for solving the algebraic system only during fault on and
356c4762a1bSJed Brown off times. It computes the entire F and then zeros out the part corresponding to
357c4762a1bSJed Brown differential equations
358c4762a1bSJed Brown F = [0;g(y)];
359c4762a1bSJed Brown */
AlgFunction(SNES snes,Vec X,Vec F,PetscCtx ctx)360*2a8381b2SBarry Smith PetscErrorCode AlgFunction(SNES snes, Vec X, Vec F, PetscCtx ctx)
361d71ae5a4SJacob Faibussowitsch {
362c4762a1bSJed Brown Userctx *user = (Userctx *)ctx;
363c4762a1bSJed Brown PetscScalar *f;
364c4762a1bSJed Brown PetscInt i;
365c4762a1bSJed Brown
366c4762a1bSJed Brown PetscFunctionBegin;
3679566063dSJacob Faibussowitsch PetscCall(ResidualFunction(snes, X, F, user));
3689566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f));
369c4762a1bSJed Brown for (i = 0; i < ngen; i++) {
370c4762a1bSJed Brown f[9 * i] = 0;
371c4762a1bSJed Brown f[9 * i + 1] = 0;
372c4762a1bSJed Brown f[9 * i + 2] = 0;
373c4762a1bSJed Brown f[9 * i + 3] = 0;
374c4762a1bSJed Brown f[9 * i + 6] = 0;
375c4762a1bSJed Brown f[9 * i + 7] = 0;
376c4762a1bSJed Brown f[9 * i + 8] = 0;
377c4762a1bSJed Brown }
3789566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f));
3793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
380c4762a1bSJed Brown }
381c4762a1bSJed Brown
PreallocateJacobian(Mat J,Userctx * user)382d71ae5a4SJacob Faibussowitsch PetscErrorCode PreallocateJacobian(Mat J, Userctx *user)
383d71ae5a4SJacob Faibussowitsch {
384c4762a1bSJed Brown PetscInt *d_nnz;
385c4762a1bSJed Brown PetscInt i, idx = 0, start = 0;
386c4762a1bSJed Brown PetscInt ncols;
387c4762a1bSJed Brown
388c4762a1bSJed Brown PetscFunctionBegin;
3899566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->neqs_pgrid, &d_nnz));
390c4762a1bSJed Brown for (i = 0; i < user->neqs_pgrid; i++) d_nnz[i] = 0;
391c4762a1bSJed Brown /* Generator subsystem */
392c4762a1bSJed Brown for (i = 0; i < ngen; i++) {
393c4762a1bSJed Brown d_nnz[idx] += 3;
394c4762a1bSJed Brown d_nnz[idx + 1] += 2;
395c4762a1bSJed Brown d_nnz[idx + 2] += 2;
396c4762a1bSJed Brown d_nnz[idx + 3] += 5;
397c4762a1bSJed Brown d_nnz[idx + 4] += 6;
398c4762a1bSJed Brown d_nnz[idx + 5] += 6;
399c4762a1bSJed Brown
400c4762a1bSJed Brown d_nnz[user->neqs_gen + 2 * gbus[i]] += 3;
401c4762a1bSJed Brown d_nnz[user->neqs_gen + 2 * gbus[i] + 1] += 3;
402c4762a1bSJed Brown
403c4762a1bSJed Brown d_nnz[idx + 6] += 2;
404c4762a1bSJed Brown d_nnz[idx + 7] += 2;
405c4762a1bSJed Brown d_nnz[idx + 8] += 5;
406c4762a1bSJed Brown
407c4762a1bSJed Brown idx = idx + 9;
408c4762a1bSJed Brown }
409c4762a1bSJed Brown
410c4762a1bSJed Brown start = user->neqs_gen;
411c4762a1bSJed Brown
412c4762a1bSJed Brown for (i = 0; i < nbus; i++) {
4139566063dSJacob Faibussowitsch PetscCall(MatGetRow(user->Ybus, 2 * i, &ncols, NULL, NULL));
414c4762a1bSJed Brown d_nnz[start + 2 * i] += ncols;
415c4762a1bSJed Brown d_nnz[start + 2 * i + 1] += ncols;
4169566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(user->Ybus, 2 * i, &ncols, NULL, NULL));
417c4762a1bSJed Brown }
418c4762a1bSJed Brown
4199566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(J, 0, d_nnz));
420c4762a1bSJed Brown
4219566063dSJacob Faibussowitsch PetscCall(PetscFree(d_nnz));
4223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
423c4762a1bSJed Brown }
424c4762a1bSJed Brown
425c4762a1bSJed Brown /*
426c4762a1bSJed Brown J = [-df_dx, -df_dy
427c4762a1bSJed Brown dg_dx, dg_dy]
428c4762a1bSJed Brown */
ResidualJacobian(SNES snes,Vec X,Mat J,Mat B,PetscCtx ctx)429*2a8381b2SBarry Smith PetscErrorCode ResidualJacobian(SNES snes, Vec X, Mat J, Mat B, PetscCtx ctx)
430d71ae5a4SJacob Faibussowitsch {
431c4762a1bSJed Brown Userctx *user = (Userctx *)ctx;
432c4762a1bSJed Brown Vec Xgen, Xnet;
433c4762a1bSJed Brown PetscScalar *xgen, *xnet;
434c4762a1bSJed Brown PetscInt i, idx = 0;
435c4762a1bSJed Brown PetscScalar Vr, Vi, Vm, Vm2;
436c4762a1bSJed Brown PetscScalar Eqp, Edp, delta; /* Generator variables */
437c4762a1bSJed Brown PetscScalar Efd; /* Exciter variables */
438c4762a1bSJed Brown PetscScalar Id, Iq; /* Generator dq axis currents */
439c4762a1bSJed Brown PetscScalar Vd, Vq;
440c4762a1bSJed Brown PetscScalar val[10];
441c4762a1bSJed Brown PetscInt row[2], col[10];
442c4762a1bSJed Brown PetscInt net_start = user->neqs_gen;
443c4762a1bSJed Brown PetscScalar Zdq_inv[4], det;
444c4762a1bSJed Brown PetscScalar dVd_dVr, dVd_dVi, dVq_dVr, dVq_dVi, dVd_ddelta, dVq_ddelta;
445c4762a1bSJed Brown PetscScalar dIGr_ddelta, dIGi_ddelta, dIGr_dId, dIGr_dIq, dIGi_dId, dIGi_dIq;
446c4762a1bSJed Brown PetscScalar dSE_dEfd;
447c4762a1bSJed Brown PetscScalar dVm_dVd, dVm_dVq, dVm_dVr, dVm_dVi;
448c4762a1bSJed Brown PetscInt ncols;
449c4762a1bSJed Brown const PetscInt *cols;
450c4762a1bSJed Brown const PetscScalar *yvals;
451c4762a1bSJed Brown PetscInt k;
452c4762a1bSJed Brown PetscScalar PD, QD, Vm0, *v0, Vm4;
453c4762a1bSJed Brown PetscScalar dPD_dVr, dPD_dVi, dQD_dVr, dQD_dVi;
454c4762a1bSJed Brown PetscScalar dIDr_dVr, dIDr_dVi, dIDi_dVr, dIDi_dVi;
455c4762a1bSJed Brown
456c4762a1bSJed Brown PetscFunctionBegin;
4579566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(B));
4589566063dSJacob Faibussowitsch PetscCall(DMCompositeGetLocalVectors(user->dmpgrid, &Xgen, &Xnet));
4599566063dSJacob Faibussowitsch PetscCall(DMCompositeScatter(user->dmpgrid, X, Xgen, Xnet));
460c4762a1bSJed Brown
4619566063dSJacob Faibussowitsch PetscCall(VecGetArray(Xgen, &xgen));
4629566063dSJacob Faibussowitsch PetscCall(VecGetArray(Xnet, &xnet));
463c4762a1bSJed Brown
464c4762a1bSJed Brown /* Generator subsystem */
465c4762a1bSJed Brown for (i = 0; i < ngen; i++) {
466c4762a1bSJed Brown Eqp = xgen[idx];
467c4762a1bSJed Brown Edp = xgen[idx + 1];
468c4762a1bSJed Brown delta = xgen[idx + 2];
469c4762a1bSJed Brown Id = xgen[idx + 4];
470c4762a1bSJed Brown Iq = xgen[idx + 5];
471c4762a1bSJed Brown Efd = xgen[idx + 6];
472c4762a1bSJed Brown
473c4762a1bSJed Brown /* fgen[idx] = (Eqp + (Xd[i] - Xdp[i])*Id - Efd)/Td0p[i]; */
474c4762a1bSJed Brown row[0] = idx;
4759371c9d4SSatish Balay col[0] = idx;
4769371c9d4SSatish Balay col[1] = idx + 4;
4779371c9d4SSatish Balay col[2] = idx + 6;
4789371c9d4SSatish Balay val[0] = 1 / Td0p[i];
4799371c9d4SSatish Balay val[1] = (Xd[i] - Xdp[i]) / Td0p[i];
4809371c9d4SSatish Balay val[2] = -1 / Td0p[i];
481c4762a1bSJed Brown
4829566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, row, 3, col, val, INSERT_VALUES));
483c4762a1bSJed Brown
484c4762a1bSJed Brown /* fgen[idx+1] = (Edp - (Xq[i] - Xqp[i])*Iq)/Tq0p[i]; */
485c4762a1bSJed Brown row[0] = idx + 1;
4869371c9d4SSatish Balay col[0] = idx + 1;
4879371c9d4SSatish Balay col[1] = idx + 5;
4889371c9d4SSatish Balay val[0] = 1 / Tq0p[i];
4899371c9d4SSatish Balay val[1] = -(Xq[i] - Xqp[i]) / Tq0p[i];
4909566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, row, 2, col, val, INSERT_VALUES));
491c4762a1bSJed Brown
492c4762a1bSJed Brown /* fgen[idx+2] = - w + w_s; */
493c4762a1bSJed Brown row[0] = idx + 2;
4949371c9d4SSatish Balay col[0] = idx + 2;
4959371c9d4SSatish Balay col[1] = idx + 3;
4969371c9d4SSatish Balay val[0] = 0;
4979371c9d4SSatish Balay val[1] = -1;
4989566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, row, 2, col, val, INSERT_VALUES));
499c4762a1bSJed Brown
500c4762a1bSJed Brown /* fgen[idx+3] = (-TM[i] + Edp*Id + Eqp*Iq + (Xqp[i] - Xdp[i])*Id*Iq + D[i]*(w - w_s))/M[i]; */
501c4762a1bSJed Brown row[0] = idx + 3;
5029371c9d4SSatish Balay col[0] = idx;
5039371c9d4SSatish Balay col[1] = idx + 1;
5049371c9d4SSatish Balay col[2] = idx + 3;
5059371c9d4SSatish Balay col[3] = idx + 4;
5069371c9d4SSatish Balay col[4] = idx + 5;
5079371c9d4SSatish Balay val[0] = Iq / M[i];
5089371c9d4SSatish Balay val[1] = Id / M[i];
5099371c9d4SSatish Balay val[2] = D[i] / M[i];
5109371c9d4SSatish Balay val[3] = (Edp + (Xqp[i] - Xdp[i]) * Iq) / M[i];
5119371c9d4SSatish Balay val[4] = (Eqp + (Xqp[i] - Xdp[i]) * Id) / M[i];
5129566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, row, 5, col, val, INSERT_VALUES));
513c4762a1bSJed Brown
514c4762a1bSJed Brown Vr = xnet[2 * gbus[i]]; /* Real part of generator terminal voltage */
515c4762a1bSJed Brown Vi = xnet[2 * gbus[i] + 1]; /* Imaginary part of the generator terminal voltage */
5169566063dSJacob Faibussowitsch PetscCall(ri2dq(Vr, Vi, delta, &Vd, &Vq));
517c4762a1bSJed Brown
518c4762a1bSJed Brown det = Rs[i] * Rs[i] + Xdp[i] * Xqp[i];
519c4762a1bSJed Brown
520c4762a1bSJed Brown Zdq_inv[0] = Rs[i] / det;
521c4762a1bSJed Brown Zdq_inv[1] = Xqp[i] / det;
522c4762a1bSJed Brown Zdq_inv[2] = -Xdp[i] / det;
523c4762a1bSJed Brown Zdq_inv[3] = Rs[i] / det;
524c4762a1bSJed Brown
5259371c9d4SSatish Balay dVd_dVr = PetscSinScalar(delta);
5269371c9d4SSatish Balay dVd_dVi = -PetscCosScalar(delta);
5279371c9d4SSatish Balay dVq_dVr = PetscCosScalar(delta);
5289371c9d4SSatish Balay dVq_dVi = PetscSinScalar(delta);
529c4762a1bSJed Brown dVd_ddelta = Vr * PetscCosScalar(delta) + Vi * PetscSinScalar(delta);
530c4762a1bSJed Brown dVq_ddelta = -Vr * PetscSinScalar(delta) + Vi * PetscCosScalar(delta);
531c4762a1bSJed Brown
532c4762a1bSJed Brown /* fgen[idx+4] = Zdq_inv[0]*(-Edp + Vd) + Zdq_inv[1]*(-Eqp + Vq) + Id; */
533c4762a1bSJed Brown row[0] = idx + 4;
5349371c9d4SSatish Balay col[0] = idx;
5359371c9d4SSatish Balay col[1] = idx + 1;
5369371c9d4SSatish Balay col[2] = idx + 2;
5379371c9d4SSatish Balay val[0] = -Zdq_inv[1];
5389371c9d4SSatish Balay val[1] = -Zdq_inv[0];
5399371c9d4SSatish Balay val[2] = Zdq_inv[0] * dVd_ddelta + Zdq_inv[1] * dVq_ddelta;
5409371c9d4SSatish Balay col[3] = idx + 4;
5419371c9d4SSatish Balay col[4] = net_start + 2 * gbus[i];
5429371c9d4SSatish Balay col[5] = net_start + 2 * gbus[i] + 1;
5439371c9d4SSatish Balay val[3] = 1;
5449371c9d4SSatish Balay val[4] = Zdq_inv[0] * dVd_dVr + Zdq_inv[1] * dVq_dVr;
5459371c9d4SSatish Balay val[5] = Zdq_inv[0] * dVd_dVi + Zdq_inv[1] * dVq_dVi;
5469566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, row, 6, col, val, INSERT_VALUES));
547c4762a1bSJed Brown
548c4762a1bSJed Brown /* fgen[idx+5] = Zdq_inv[2]*(-Edp + Vd) + Zdq_inv[3]*(-Eqp + Vq) + Iq; */
549c4762a1bSJed Brown row[0] = idx + 5;
5509371c9d4SSatish Balay col[0] = idx;
5519371c9d4SSatish Balay col[1] = idx + 1;
5529371c9d4SSatish Balay col[2] = idx + 2;
5539371c9d4SSatish Balay val[0] = -Zdq_inv[3];
5549371c9d4SSatish Balay val[1] = -Zdq_inv[2];
5559371c9d4SSatish Balay val[2] = Zdq_inv[2] * dVd_ddelta + Zdq_inv[3] * dVq_ddelta;
5569371c9d4SSatish Balay col[3] = idx + 5;
5579371c9d4SSatish Balay col[4] = net_start + 2 * gbus[i];
5589371c9d4SSatish Balay col[5] = net_start + 2 * gbus[i] + 1;
5599371c9d4SSatish Balay val[3] = 1;
5609371c9d4SSatish Balay val[4] = Zdq_inv[2] * dVd_dVr + Zdq_inv[3] * dVq_dVr;
5619371c9d4SSatish Balay val[5] = Zdq_inv[2] * dVd_dVi + Zdq_inv[3] * dVq_dVi;
5629566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, row, 6, col, val, INSERT_VALUES));
563c4762a1bSJed Brown
564c4762a1bSJed Brown dIGr_ddelta = Id * PetscCosScalar(delta) - Iq * PetscSinScalar(delta);
565c4762a1bSJed Brown dIGi_ddelta = Id * PetscSinScalar(delta) + Iq * PetscCosScalar(delta);
5669371c9d4SSatish Balay dIGr_dId = PetscSinScalar(delta);
5679371c9d4SSatish Balay dIGr_dIq = PetscCosScalar(delta);
5689371c9d4SSatish Balay dIGi_dId = -PetscCosScalar(delta);
5699371c9d4SSatish Balay dIGi_dIq = PetscSinScalar(delta);
570c4762a1bSJed Brown
571c4762a1bSJed Brown /* fnet[2*gbus[i]] -= IGi; */
572c4762a1bSJed Brown row[0] = net_start + 2 * gbus[i];
5739371c9d4SSatish Balay col[0] = idx + 2;
5749371c9d4SSatish Balay col[1] = idx + 4;
5759371c9d4SSatish Balay col[2] = idx + 5;
5769371c9d4SSatish Balay val[0] = -dIGi_ddelta;
5779371c9d4SSatish Balay val[1] = -dIGi_dId;
5789371c9d4SSatish Balay val[2] = -dIGi_dIq;
5799566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, row, 3, col, val, INSERT_VALUES));
580c4762a1bSJed Brown
581c4762a1bSJed Brown /* fnet[2*gbus[i]+1] -= IGr; */
582c4762a1bSJed Brown row[0] = net_start + 2 * gbus[i] + 1;
5839371c9d4SSatish Balay col[0] = idx + 2;
5849371c9d4SSatish Balay col[1] = idx + 4;
5859371c9d4SSatish Balay col[2] = idx + 5;
5869371c9d4SSatish Balay val[0] = -dIGr_ddelta;
5879371c9d4SSatish Balay val[1] = -dIGr_dId;
5889371c9d4SSatish Balay val[2] = -dIGr_dIq;
5899566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, row, 3, col, val, INSERT_VALUES));
590c4762a1bSJed Brown
591c4762a1bSJed Brown Vm = PetscSqrtScalar(Vd * Vd + Vq * Vq);
592c4762a1bSJed Brown
593c4762a1bSJed Brown /* fgen[idx+6] = (KE[i]*Efd + SE - VR)/TE[i]; */
594c4762a1bSJed Brown /* SE = k1[i]*PetscExpScalar(k2[i]*Efd); */
595c4762a1bSJed Brown
596c4762a1bSJed Brown dSE_dEfd = k1[i] * k2[i] * PetscExpScalar(k2[i] * Efd);
597c4762a1bSJed Brown
598c4762a1bSJed Brown row[0] = idx + 6;
5999371c9d4SSatish Balay col[0] = idx + 6;
6009371c9d4SSatish Balay col[1] = idx + 8;
6019371c9d4SSatish Balay val[0] = (KE[i] + dSE_dEfd) / TE[i];
6029371c9d4SSatish Balay val[1] = -1 / TE[i];
6039566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, row, 2, col, val, INSERT_VALUES));
604c4762a1bSJed Brown
605c4762a1bSJed Brown /* Exciter differential equations */
606c4762a1bSJed Brown
607c4762a1bSJed Brown /* fgen[idx+7] = (RF - KF[i]*Efd/TF[i])/TF[i]; */
608c4762a1bSJed Brown row[0] = idx + 7;
6099371c9d4SSatish Balay col[0] = idx + 6;
6109371c9d4SSatish Balay col[1] = idx + 7;
6119371c9d4SSatish Balay val[0] = (-KF[i] / TF[i]) / TF[i];
6129371c9d4SSatish Balay val[1] = 1 / TF[i];
6139566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, row, 2, col, val, INSERT_VALUES));
614c4762a1bSJed Brown
615c4762a1bSJed Brown /* fgen[idx+8] = (VR - KA[i]*RF + KA[i]*KF[i]*Efd/TF[i] - KA[i]*(Vref[i] - Vm))/TA[i]; */
616c4762a1bSJed Brown /* Vm = (Vd^2 + Vq^2)^0.5; */
6179371c9d4SSatish Balay dVm_dVd = Vd / Vm;
6189371c9d4SSatish Balay dVm_dVq = Vq / Vm;
619c4762a1bSJed Brown dVm_dVr = dVm_dVd * dVd_dVr + dVm_dVq * dVq_dVr;
620c4762a1bSJed Brown dVm_dVi = dVm_dVd * dVd_dVi + dVm_dVq * dVq_dVi;
621c4762a1bSJed Brown row[0] = idx + 8;
6229371c9d4SSatish Balay col[0] = idx + 6;
6239371c9d4SSatish Balay col[1] = idx + 7;
6249371c9d4SSatish Balay col[2] = idx + 8;
6259371c9d4SSatish Balay val[0] = (KA[i] * KF[i] / TF[i]) / TA[i];
6269371c9d4SSatish Balay val[1] = -KA[i] / TA[i];
6279371c9d4SSatish Balay val[2] = 1 / TA[i];
6289371c9d4SSatish Balay col[3] = net_start + 2 * gbus[i];
6299371c9d4SSatish Balay col[4] = net_start + 2 * gbus[i] + 1;
6309371c9d4SSatish Balay val[3] = KA[i] * dVm_dVr / TA[i];
6319371c9d4SSatish Balay val[4] = KA[i] * dVm_dVi / TA[i];
6329566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, row, 5, col, val, INSERT_VALUES));
633c4762a1bSJed Brown idx = idx + 9;
634c4762a1bSJed Brown }
635c4762a1bSJed Brown
636c4762a1bSJed Brown for (i = 0; i < nbus; i++) {
6379566063dSJacob Faibussowitsch PetscCall(MatGetRow(user->Ybus, 2 * i, &ncols, &cols, &yvals));
638c4762a1bSJed Brown row[0] = net_start + 2 * i;
639c4762a1bSJed Brown for (k = 0; k < ncols; k++) {
640c4762a1bSJed Brown col[k] = net_start + cols[k];
641c4762a1bSJed Brown val[k] = yvals[k];
642c4762a1bSJed Brown }
6439566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, row, ncols, col, val, INSERT_VALUES));
6449566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(user->Ybus, 2 * i, &ncols, &cols, &yvals));
645c4762a1bSJed Brown
6469566063dSJacob Faibussowitsch PetscCall(MatGetRow(user->Ybus, 2 * i + 1, &ncols, &cols, &yvals));
647c4762a1bSJed Brown row[0] = net_start + 2 * i + 1;
648c4762a1bSJed Brown for (k = 0; k < ncols; k++) {
649c4762a1bSJed Brown col[k] = net_start + cols[k];
650c4762a1bSJed Brown val[k] = yvals[k];
651c4762a1bSJed Brown }
6529566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, row, ncols, col, val, INSERT_VALUES));
6539566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(user->Ybus, 2 * i + 1, &ncols, &cols, &yvals));
654c4762a1bSJed Brown }
655c4762a1bSJed Brown
6569566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J, MAT_FLUSH_ASSEMBLY));
6579566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J, MAT_FLUSH_ASSEMBLY));
658c4762a1bSJed Brown
6599566063dSJacob Faibussowitsch PetscCall(VecGetArray(user->V0, &v0));
660c4762a1bSJed Brown for (i = 0; i < nload; i++) {
661c4762a1bSJed Brown Vr = xnet[2 * lbus[i]]; /* Real part of load bus voltage */
662c4762a1bSJed Brown Vi = xnet[2 * lbus[i] + 1]; /* Imaginary part of the load bus voltage */
6639371c9d4SSatish Balay Vm = PetscSqrtScalar(Vr * Vr + Vi * Vi);
6649371c9d4SSatish Balay Vm2 = Vm * Vm;
6659371c9d4SSatish Balay Vm4 = Vm2 * Vm2;
666c4762a1bSJed Brown Vm0 = PetscSqrtScalar(v0[2 * lbus[i]] * v0[2 * lbus[i]] + v0[2 * lbus[i] + 1] * v0[2 * lbus[i] + 1]);
667c4762a1bSJed Brown PD = QD = 0.0;
668c4762a1bSJed Brown dPD_dVr = dPD_dVi = dQD_dVr = dQD_dVi = 0.0;
669c4762a1bSJed Brown for (k = 0; k < ld_nsegsp[i]; k++) {
67057508eceSPierre Jolivet PD += ld_alphap[k] * PD0[i] * PetscPowScalar(Vm / Vm0, ld_betap[k]);
67157508eceSPierre Jolivet dPD_dVr += ld_alphap[k] * ld_betap[k] * PD0[i] * PetscPowScalar(1 / Vm0, ld_betap[k]) * Vr * PetscPowScalar(Vm, ld_betap[k] - 2);
67257508eceSPierre Jolivet dPD_dVi += ld_alphap[k] * ld_betap[k] * PD0[i] * PetscPowScalar(1 / Vm0, ld_betap[k]) * Vi * PetscPowScalar(Vm, ld_betap[k] - 2);
673c4762a1bSJed Brown }
674c4762a1bSJed Brown for (k = 0; k < ld_nsegsq[i]; k++) {
67557508eceSPierre Jolivet QD += ld_alphaq[k] * QD0[i] * PetscPowScalar(Vm / Vm0, ld_betaq[k]);
67657508eceSPierre Jolivet dQD_dVr += ld_alphaq[k] * ld_betaq[k] * QD0[i] * PetscPowScalar(1 / Vm0, ld_betaq[k]) * Vr * PetscPowScalar(Vm, ld_betaq[k] - 2);
67757508eceSPierre Jolivet dQD_dVi += ld_alphaq[k] * ld_betaq[k] * QD0[i] * PetscPowScalar(1 / Vm0, ld_betaq[k]) * Vi * PetscPowScalar(Vm, ld_betaq[k] - 2);
678c4762a1bSJed Brown }
679c4762a1bSJed Brown
680c4762a1bSJed Brown /* IDr = (PD*Vr + QD*Vi)/Vm2; */
681c4762a1bSJed Brown /* IDi = (-QD*Vr + PD*Vi)/Vm2; */
682c4762a1bSJed Brown
683c4762a1bSJed Brown dIDr_dVr = (dPD_dVr * Vr + dQD_dVr * Vi + PD) / Vm2 - ((PD * Vr + QD * Vi) * 2 * Vr) / Vm4;
684c4762a1bSJed Brown dIDr_dVi = (dPD_dVi * Vr + dQD_dVi * Vi + QD) / Vm2 - ((PD * Vr + QD * Vi) * 2 * Vi) / Vm4;
685c4762a1bSJed Brown
686c4762a1bSJed Brown dIDi_dVr = (-dQD_dVr * Vr + dPD_dVr * Vi - QD) / Vm2 - ((-QD * Vr + PD * Vi) * 2 * Vr) / Vm4;
687c4762a1bSJed Brown dIDi_dVi = (-dQD_dVi * Vr + dPD_dVi * Vi + PD) / Vm2 - ((-QD * Vr + PD * Vi) * 2 * Vi) / Vm4;
688c4762a1bSJed Brown
689c4762a1bSJed Brown /* fnet[2*lbus[i]] += IDi; */
690c4762a1bSJed Brown row[0] = net_start + 2 * lbus[i];
6919371c9d4SSatish Balay col[0] = net_start + 2 * lbus[i];
6929371c9d4SSatish Balay col[1] = net_start + 2 * lbus[i] + 1;
6939371c9d4SSatish Balay val[0] = dIDi_dVr;
6949371c9d4SSatish Balay val[1] = dIDi_dVi;
6959566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, row, 2, col, val, ADD_VALUES));
696c4762a1bSJed Brown /* fnet[2*lbus[i]+1] += IDr; */
697c4762a1bSJed Brown row[0] = net_start + 2 * lbus[i] + 1;
6989371c9d4SSatish Balay col[0] = net_start + 2 * lbus[i];
6999371c9d4SSatish Balay col[1] = net_start + 2 * lbus[i] + 1;
7009371c9d4SSatish Balay val[0] = dIDr_dVr;
7019371c9d4SSatish Balay val[1] = dIDr_dVi;
7029566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, row, 2, col, val, ADD_VALUES));
703c4762a1bSJed Brown }
7049566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(user->V0, &v0));
705c4762a1bSJed Brown
7069566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Xgen, &xgen));
7079566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Xnet, &xnet));
708c4762a1bSJed Brown
7099566063dSJacob Faibussowitsch PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid, &Xgen, &Xnet));
710c4762a1bSJed Brown
7119566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
7129566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
7133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
714c4762a1bSJed Brown }
715c4762a1bSJed Brown
716c4762a1bSJed Brown /*
717c4762a1bSJed Brown J = [I, 0
718c4762a1bSJed Brown dg_dx, dg_dy]
719c4762a1bSJed Brown */
AlgJacobian(SNES snes,Vec X,Mat A,Mat B,PetscCtx ctx)720*2a8381b2SBarry Smith PetscErrorCode AlgJacobian(SNES snes, Vec X, Mat A, Mat B, PetscCtx ctx)
721d71ae5a4SJacob Faibussowitsch {
722c4762a1bSJed Brown Userctx *user = (Userctx *)ctx;
723c4762a1bSJed Brown
724c4762a1bSJed Brown PetscFunctionBegin;
7259566063dSJacob Faibussowitsch PetscCall(ResidualJacobian(snes, X, A, B, ctx));
7269566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_KEEP_NONZERO_PATTERN, PETSC_TRUE));
7279566063dSJacob Faibussowitsch PetscCall(MatZeroRowsIS(A, user->is_diff, 1.0, NULL, NULL));
7283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
729c4762a1bSJed Brown }
730c4762a1bSJed Brown
731c4762a1bSJed Brown /*
732c4762a1bSJed Brown J = [a*I-df_dx, -df_dy
733c4762a1bSJed Brown dg_dx, dg_dy]
734c4762a1bSJed Brown */
735c4762a1bSJed Brown
IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,Userctx * user)736d71ae5a4SJacob Faibussowitsch PetscErrorCode IJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat A, Mat B, Userctx *user)
737d71ae5a4SJacob Faibussowitsch {
738c4762a1bSJed Brown SNES snes;
739c4762a1bSJed Brown PetscScalar atmp = (PetscScalar)a;
740c4762a1bSJed Brown PetscInt i, row;
741c4762a1bSJed Brown
742c4762a1bSJed Brown PetscFunctionBegin;
743c4762a1bSJed Brown user->t = t;
744c4762a1bSJed Brown
7459566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes));
7469566063dSJacob Faibussowitsch PetscCall(ResidualJacobian(snes, X, A, B, user));
747c4762a1bSJed Brown for (i = 0; i < ngen; i++) {
748c4762a1bSJed Brown row = 9 * i;
7499566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES));
750c4762a1bSJed Brown row = 9 * i + 1;
7519566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES));
752c4762a1bSJed Brown row = 9 * i + 2;
7539566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES));
754c4762a1bSJed Brown row = 9 * i + 3;
7559566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES));
756c4762a1bSJed Brown row = 9 * i + 6;
7579566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES));
758c4762a1bSJed Brown row = 9 * i + 7;
7599566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES));
760c4762a1bSJed Brown row = 9 * i + 8;
7619566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES));
762c4762a1bSJed Brown }
7639566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
7649566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
7653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
766c4762a1bSJed Brown }
767c4762a1bSJed Brown
CostIntegrand(TS ts,PetscReal t,Vec U,Vec R,Userctx * user)768d71ae5a4SJacob Faibussowitsch static PetscErrorCode CostIntegrand(TS ts, PetscReal t, Vec U, Vec R, Userctx *user)
769d71ae5a4SJacob Faibussowitsch {
770c4762a1bSJed Brown PetscScalar *r;
771c4762a1bSJed Brown const PetscScalar *u;
772c4762a1bSJed Brown PetscInt idx;
773c4762a1bSJed Brown Vec Xgen, Xnet;
774c4762a1bSJed Brown PetscScalar *xgen;
775c4762a1bSJed Brown PetscInt i;
776c4762a1bSJed Brown
777c4762a1bSJed Brown PetscFunctionBegin;
7789566063dSJacob Faibussowitsch PetscCall(DMCompositeGetLocalVectors(user->dmpgrid, &Xgen, &Xnet));
7799566063dSJacob Faibussowitsch PetscCall(DMCompositeScatter(user->dmpgrid, U, Xgen, Xnet));
780c4762a1bSJed Brown
7819566063dSJacob Faibussowitsch PetscCall(VecGetArray(Xgen, &xgen));
782c4762a1bSJed Brown
7839566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u));
7849566063dSJacob Faibussowitsch PetscCall(VecGetArray(R, &r));
785c4762a1bSJed Brown r[0] = 0.;
786c4762a1bSJed Brown
787c4762a1bSJed Brown idx = 0;
788c4762a1bSJed Brown for (i = 0; i < ngen; i++) {
789c4762a1bSJed 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);
790c4762a1bSJed Brown idx += 9;
791c4762a1bSJed Brown }
7929566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(R, &r));
7939566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u));
7949566063dSJacob Faibussowitsch PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid, &Xgen, &Xnet));
7953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
796c4762a1bSJed Brown }
797c4762a1bSJed Brown
MonitorUpdateQ(TS ts,PetscInt stepnum,PetscReal time,Vec X,PetscCtx ctx0)798*2a8381b2SBarry Smith static PetscErrorCode MonitorUpdateQ(TS ts, PetscInt stepnum, PetscReal time, Vec X, PetscCtx ctx0)
799d71ae5a4SJacob Faibussowitsch {
800c4762a1bSJed Brown Vec C, *Y;
801c4762a1bSJed Brown PetscInt Nr;
802c4762a1bSJed Brown PetscReal h, theta;
803c4762a1bSJed Brown Userctx *ctx = (Userctx *)ctx0;
804c4762a1bSJed Brown
805c4762a1bSJed Brown PetscFunctionBegin;
806c4762a1bSJed Brown theta = 0.5;
8079566063dSJacob Faibussowitsch PetscCall(TSGetStages(ts, &Nr, &Y));
8089566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts, &h));
8099566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ctx->vec_q, &C));
810c4762a1bSJed Brown /* compute integrals */
811c4762a1bSJed Brown if (stepnum > 0) {
8129566063dSJacob Faibussowitsch PetscCall(CostIntegrand(ts, time, X, C, ctx));
8139566063dSJacob Faibussowitsch PetscCall(VecAXPY(ctx->vec_q, h * theta, C));
8149566063dSJacob Faibussowitsch PetscCall(CostIntegrand(ts, time + h * theta, Y[0], C, ctx));
8159566063dSJacob Faibussowitsch PetscCall(VecAXPY(ctx->vec_q, h * (1 - theta), C));
816c4762a1bSJed Brown }
8179566063dSJacob Faibussowitsch PetscCall(VecDestroy(&C));
8183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
819c4762a1bSJed Brown }
820c4762a1bSJed Brown
main(int argc,char ** argv)821d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
822d71ae5a4SJacob Faibussowitsch {
823c4762a1bSJed Brown Userctx user;
824c4762a1bSJed Brown Vec p;
825c4762a1bSJed Brown PetscScalar *x_ptr;
826c4762a1bSJed Brown PetscMPIInt size;
827c4762a1bSJed Brown PetscInt i;
828c4762a1bSJed Brown KSP ksp;
829c4762a1bSJed Brown PC pc;
830c4762a1bSJed Brown PetscInt *idx2;
831c4762a1bSJed Brown Tao tao;
832c4762a1bSJed Brown Vec lowerb, upperb;
833c4762a1bSJed Brown
834c4762a1bSJed Brown PetscFunctionBeginUser;
835327415f7SBarry Smith PetscFunctionBeginUser;
8369566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, "petscoptions", help));
8379566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
8383c633725SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs");
839c4762a1bSJed Brown
8409566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_WORLD, 1, &user.vec_q));
841c4762a1bSJed Brown
842c4762a1bSJed Brown user.neqs_gen = 9 * ngen; /* # eqs. for generator subsystem */
843c4762a1bSJed Brown user.neqs_net = 2 * nbus; /* # eqs. for network subsystem */
844c4762a1bSJed Brown user.neqs_pgrid = user.neqs_gen + user.neqs_net;
845c4762a1bSJed Brown
846c4762a1bSJed Brown /* Create indices for differential and algebraic equations */
8479566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(7 * ngen, &idx2));
848c4762a1bSJed Brown for (i = 0; i < ngen; i++) {
8499371c9d4SSatish Balay idx2[7 * i] = 9 * i;
8509371c9d4SSatish Balay idx2[7 * i + 1] = 9 * i + 1;
8519371c9d4SSatish Balay idx2[7 * i + 2] = 9 * i + 2;
8529371c9d4SSatish Balay idx2[7 * i + 3] = 9 * i + 3;
8539371c9d4SSatish Balay idx2[7 * i + 4] = 9 * i + 6;
8549371c9d4SSatish Balay idx2[7 * i + 5] = 9 * i + 7;
8559371c9d4SSatish Balay idx2[7 * i + 6] = 9 * i + 8;
856c4762a1bSJed Brown }
8579566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, 7 * ngen, idx2, PETSC_COPY_VALUES, &user.is_diff));
8589566063dSJacob Faibussowitsch PetscCall(ISComplement(user.is_diff, 0, user.neqs_pgrid, &user.is_alg));
8599566063dSJacob Faibussowitsch PetscCall(PetscFree(idx2));
860c4762a1bSJed Brown
861c4762a1bSJed Brown /* Set run time options */
862d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Transient stability fault options", "");
863c4762a1bSJed Brown {
864c4762a1bSJed Brown user.tfaulton = 1.0;
865c4762a1bSJed Brown user.tfaultoff = 1.2;
866c4762a1bSJed Brown user.Rfault = 0.0001;
867c4762a1bSJed Brown user.faultbus = 8;
8689566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tfaulton", "", "", user.tfaulton, &user.tfaulton, NULL));
8699566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tfaultoff", "", "", user.tfaultoff, &user.tfaultoff, NULL));
8709566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-faultbus", "", "", user.faultbus, &user.faultbus, NULL));
871c4762a1bSJed Brown user.t0 = 0.0;
872c4762a1bSJed Brown user.tmax = 1.5;
8739566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-t0", "", "", user.t0, &user.t0, NULL));
8749566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tmax", "", "", user.tmax, &user.tmax, NULL));
875c4762a1bSJed Brown user.freq_u = 61.0;
876c4762a1bSJed Brown user.freq_l = 59.0;
877c4762a1bSJed Brown user.pow = 2;
8789566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-frequ", "", "", user.freq_u, &user.freq_u, NULL));
8799566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-freql", "", "", user.freq_l, &user.freq_l, NULL));
8809566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pow", "", "", user.pow, &user.pow, NULL));
881c4762a1bSJed Brown }
882d0609cedSBarry Smith PetscOptionsEnd();
883c4762a1bSJed Brown
884c4762a1bSJed Brown /* Create DMs for generator and network subsystems */
8859566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, user.neqs_gen, 1, 1, NULL, &user.dmgen));
8869566063dSJacob Faibussowitsch PetscCall(DMSetOptionsPrefix(user.dmgen, "dmgen_"));
8879566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(user.dmgen));
8889566063dSJacob Faibussowitsch PetscCall(DMSetUp(user.dmgen));
8899566063dSJacob Faibussowitsch PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, user.neqs_net, 1, 1, NULL, &user.dmnet));
8909566063dSJacob Faibussowitsch PetscCall(DMSetOptionsPrefix(user.dmnet, "dmnet_"));
8919566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(user.dmnet));
8929566063dSJacob Faibussowitsch PetscCall(DMSetUp(user.dmnet));
893c4762a1bSJed Brown /* Create a composite DM packer and add the two DMs */
8949566063dSJacob Faibussowitsch PetscCall(DMCompositeCreate(PETSC_COMM_WORLD, &user.dmpgrid));
8959566063dSJacob Faibussowitsch PetscCall(DMSetOptionsPrefix(user.dmpgrid, "pgrid_"));
8969566063dSJacob Faibussowitsch PetscCall(DMCompositeAddDM(user.dmpgrid, user.dmgen));
8979566063dSJacob Faibussowitsch PetscCall(DMCompositeAddDM(user.dmpgrid, user.dmnet));
898c4762a1bSJed Brown
899c4762a1bSJed Brown /* Create TAO solver and set desired solution method */
9009566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao));
9019566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao, TAOBLMVM));
902c4762a1bSJed Brown /*
903c4762a1bSJed Brown Optimization starts
904c4762a1bSJed Brown */
905c4762a1bSJed Brown /* Set initial solution guess */
9069566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_WORLD, 3, &p));
9079566063dSJacob Faibussowitsch PetscCall(VecGetArray(p, &x_ptr));
9089371c9d4SSatish Balay x_ptr[0] = PG[0];
9099371c9d4SSatish Balay x_ptr[1] = PG[1];
9109371c9d4SSatish Balay x_ptr[2] = PG[2];
9119566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(p, &x_ptr));
912c4762a1bSJed Brown
9139566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao, p));
914c4762a1bSJed Brown /* Set routine for function and gradient evaluation */
9159566063dSJacob Faibussowitsch PetscCall(TaoSetObjective(tao, FormFunction, (void *)&user));
9169566063dSJacob Faibussowitsch PetscCall(TaoSetGradient(tao, NULL, TaoDefaultComputeGradient, (void *)&user));
917c4762a1bSJed Brown
918c4762a1bSJed Brown /* Set bounds for the optimization */
9199566063dSJacob Faibussowitsch PetscCall(VecDuplicate(p, &lowerb));
9209566063dSJacob Faibussowitsch PetscCall(VecDuplicate(p, &upperb));
9219566063dSJacob Faibussowitsch PetscCall(VecGetArray(lowerb, &x_ptr));
9229371c9d4SSatish Balay x_ptr[0] = 0.5;
9239371c9d4SSatish Balay x_ptr[1] = 0.5;
9249371c9d4SSatish Balay x_ptr[2] = 0.5;
9259566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(lowerb, &x_ptr));
9269566063dSJacob Faibussowitsch PetscCall(VecGetArray(upperb, &x_ptr));
9279371c9d4SSatish Balay x_ptr[0] = 2.0;
9289371c9d4SSatish Balay x_ptr[1] = 2.0;
9299371c9d4SSatish Balay x_ptr[2] = 2.0;
9309566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(upperb, &x_ptr));
9319566063dSJacob Faibussowitsch PetscCall(TaoSetVariableBounds(tao, lowerb, upperb));
932c4762a1bSJed Brown
933c4762a1bSJed Brown /* Check for any TAO command line options */
9349566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao));
9359566063dSJacob Faibussowitsch PetscCall(TaoGetKSP(tao, &ksp));
936c4762a1bSJed Brown if (ksp) {
9379566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp, &pc));
9389566063dSJacob Faibussowitsch PetscCall(PCSetType(pc, PCNONE));
939c4762a1bSJed Brown }
940c4762a1bSJed Brown
941c4762a1bSJed Brown /* SOLVE THE APPLICATION */
9429566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao));
943c4762a1bSJed Brown
9449566063dSJacob Faibussowitsch PetscCall(VecView(p, PETSC_VIEWER_STDOUT_WORLD));
945c4762a1bSJed Brown /* Free TAO data structures */
9469566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao));
9479566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user.vec_q));
9489566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lowerb));
9499566063dSJacob Faibussowitsch PetscCall(VecDestroy(&upperb));
9509566063dSJacob Faibussowitsch PetscCall(VecDestroy(&p));
9519566063dSJacob Faibussowitsch PetscCall(DMDestroy(&user.dmgen));
9529566063dSJacob Faibussowitsch PetscCall(DMDestroy(&user.dmnet));
9539566063dSJacob Faibussowitsch PetscCall(DMDestroy(&user.dmpgrid));
9549566063dSJacob Faibussowitsch PetscCall(ISDestroy(&user.is_diff));
9559566063dSJacob Faibussowitsch PetscCall(ISDestroy(&user.is_alg));
9569566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
957b122ec5aSJacob Faibussowitsch return 0;
958c4762a1bSJed Brown }
959c4762a1bSJed Brown
960c4762a1bSJed Brown /* ------------------------------------------------------------------ */
961c4762a1bSJed Brown /*
962c4762a1bSJed Brown FormFunction - Evaluates the function and corresponding gradient.
963c4762a1bSJed Brown
964c4762a1bSJed Brown Input Parameters:
965c4762a1bSJed Brown tao - the Tao context
966c4762a1bSJed Brown X - the input vector
967a82e8c82SStefano Zampini ptr - optional user-defined context, as set by TaoSetObjectiveAndGradient()
968c4762a1bSJed Brown
969c4762a1bSJed Brown Output Parameters:
970c4762a1bSJed Brown f - the newly evaluated function
971c4762a1bSJed Brown */
FormFunction(Tao tao,Vec P,PetscReal * f,PetscCtx ctx0)972*2a8381b2SBarry Smith PetscErrorCode FormFunction(Tao tao, Vec P, PetscReal *f, PetscCtx ctx0)
973d71ae5a4SJacob Faibussowitsch {
974c4762a1bSJed Brown TS ts;
975c4762a1bSJed Brown SNES snes_alg;
976c4762a1bSJed Brown Userctx *ctx = (Userctx *)ctx0;
977c4762a1bSJed Brown Vec X;
978c4762a1bSJed Brown Mat J;
979c4762a1bSJed Brown /* sensitivity context */
980c4762a1bSJed Brown PetscScalar *x_ptr;
981c4762a1bSJed Brown PetscViewer Xview, Ybusview;
982c4762a1bSJed Brown Vec F_alg;
983c4762a1bSJed Brown Vec Xdot;
984c4762a1bSJed Brown PetscInt row_loc, col_loc;
985c4762a1bSJed Brown PetscScalar val;
986c4762a1bSJed Brown
9879566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(P, (const PetscScalar **)&x_ptr));
988c4762a1bSJed Brown PG[0] = x_ptr[0];
989c4762a1bSJed Brown PG[1] = x_ptr[1];
990c4762a1bSJed Brown PG[2] = x_ptr[2];
9919566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(P, (const PetscScalar **)&x_ptr));
992c4762a1bSJed Brown
993c4762a1bSJed Brown ctx->stepnum = 0;
994c4762a1bSJed Brown
9959566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(ctx->vec_q));
996c4762a1bSJed Brown
997c4762a1bSJed Brown /* Read initial voltage vector and Ybus */
9989566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "X.bin", FILE_MODE_READ, &Xview));
9999566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "Ybus.bin", FILE_MODE_READ, &Ybusview));
1000c4762a1bSJed Brown
10019566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &ctx->V0));
10029566063dSJacob Faibussowitsch PetscCall(VecSetSizes(ctx->V0, PETSC_DECIDE, ctx->neqs_net));
10039566063dSJacob Faibussowitsch PetscCall(VecLoad(ctx->V0, Xview));
1004c4762a1bSJed Brown
10059566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &ctx->Ybus));
10069566063dSJacob Faibussowitsch PetscCall(MatSetSizes(ctx->Ybus, PETSC_DECIDE, PETSC_DECIDE, ctx->neqs_net, ctx->neqs_net));
10079566063dSJacob Faibussowitsch PetscCall(MatSetType(ctx->Ybus, MATBAIJ));
10089566063dSJacob Faibussowitsch /* PetscCall(MatSetBlockSize(ctx->Ybus,2)); */
10099566063dSJacob Faibussowitsch PetscCall(MatLoad(ctx->Ybus, Ybusview));
1010c4762a1bSJed Brown
10119566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&Xview));
10129566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&Ybusview));
1013c4762a1bSJed Brown
10149566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(ctx->dmpgrid, &X));
1015c4762a1bSJed Brown
10169566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &J));
10179566063dSJacob Faibussowitsch PetscCall(MatSetSizes(J, PETSC_DECIDE, PETSC_DECIDE, ctx->neqs_pgrid, ctx->neqs_pgrid));
10189566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(J));
10199566063dSJacob Faibussowitsch PetscCall(PreallocateJacobian(J, ctx));
1020c4762a1bSJed Brown
1021c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1022c4762a1bSJed Brown Create timestepping solver context
1023c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
10249566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
10259566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
10269566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSCN));
10278434afd1SBarry Smith PetscCall(TSSetIFunction(ts, NULL, (TSIFunctionFn *)IFunction, ctx));
10288434afd1SBarry Smith PetscCall(TSSetIJacobian(ts, J, J, (TSIJacobianFn *)IJacobian, ctx));
10299566063dSJacob Faibussowitsch PetscCall(TSSetApplicationContext(ts, ctx));
1030c4762a1bSJed Brown
10319566063dSJacob Faibussowitsch PetscCall(TSMonitorSet(ts, MonitorUpdateQ, ctx, NULL));
1032c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1033c4762a1bSJed Brown Set initial conditions
1034c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
10359566063dSJacob Faibussowitsch PetscCall(SetInitialGuess(X, ctx));
1036c4762a1bSJed Brown
10379566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X, &F_alg));
10389566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes_alg));
10399566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes_alg, F_alg, AlgFunction, ctx));
10409566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(J));
10419566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes_alg, J, J, AlgJacobian, ctx));
10429566063dSJacob Faibussowitsch PetscCall(SNESSetOptionsPrefix(snes_alg, "alg_"));
10439566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes_alg));
1044c4762a1bSJed Brown ctx->alg_flg = PETSC_TRUE;
1045c4762a1bSJed Brown /* Solve the algebraic equations */
10469566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes_alg, NULL, X));
1047c4762a1bSJed Brown
1048c4762a1bSJed Brown /* Just to set up the Jacobian structure */
10499566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X, &Xdot));
10509566063dSJacob Faibussowitsch PetscCall(IJacobian(ts, 0.0, X, Xdot, 0.0, J, J, ctx));
10519566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Xdot));
1052c4762a1bSJed Brown
1053c4762a1bSJed Brown ctx->stepnum++;
1054c4762a1bSJed Brown
10559566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, 0.01));
10569566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, ctx->tmax));
10579566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
10589566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts));
1059c4762a1bSJed Brown
1060c4762a1bSJed Brown /* Prefault period */
1061c4762a1bSJed Brown ctx->alg_flg = PETSC_FALSE;
10629566063dSJacob Faibussowitsch PetscCall(TSSetTime(ts, 0.0));
10639566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, ctx->tfaulton));
10649566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, X));
1065c4762a1bSJed Brown
1066c4762a1bSJed Brown /* Create the nonlinear solver for solving the algebraic system */
1067c4762a1bSJed Brown /* Note that although the algebraic system needs to be solved only for
1068c4762a1bSJed Brown Idq and V, we reuse the entire system including xgen. The xgen
1069c4762a1bSJed Brown variables are held constant by setting their residuals to 0 and
1070c4762a1bSJed Brown putting a 1 on the Jacobian diagonal for xgen rows
1071c4762a1bSJed Brown */
10729566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(J));
1073c4762a1bSJed Brown
1074c4762a1bSJed Brown /* Apply disturbance - resistive fault at ctx->faultbus */
1075c4762a1bSJed Brown /* This is done by adding shunt conductance to the diagonal location
1076c4762a1bSJed Brown in the Ybus matrix */
10779371c9d4SSatish Balay row_loc = 2 * ctx->faultbus;
10789371c9d4SSatish Balay col_loc = 2 * ctx->faultbus + 1; /* Location for G */
1079c4762a1bSJed Brown val = 1 / ctx->Rfault;
10809566063dSJacob Faibussowitsch PetscCall(MatSetValues(ctx->Ybus, 1, &row_loc, 1, &col_loc, &val, ADD_VALUES));
10819371c9d4SSatish Balay row_loc = 2 * ctx->faultbus + 1;
10829371c9d4SSatish Balay col_loc = 2 * ctx->faultbus; /* Location for G */
1083c4762a1bSJed Brown val = 1 / ctx->Rfault;
10849566063dSJacob Faibussowitsch PetscCall(MatSetValues(ctx->Ybus, 1, &row_loc, 1, &col_loc, &val, ADD_VALUES));
1085c4762a1bSJed Brown
10869566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(ctx->Ybus, MAT_FINAL_ASSEMBLY));
10879566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(ctx->Ybus, MAT_FINAL_ASSEMBLY));
1088c4762a1bSJed Brown
1089c4762a1bSJed Brown ctx->alg_flg = PETSC_TRUE;
1090c4762a1bSJed Brown /* Solve the algebraic equations */
10919566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes_alg, NULL, X));
1092c4762a1bSJed Brown
1093c4762a1bSJed Brown ctx->stepnum++;
1094c4762a1bSJed Brown
1095c4762a1bSJed Brown /* Disturbance period */
1096c4762a1bSJed Brown ctx->alg_flg = PETSC_FALSE;
10979566063dSJacob Faibussowitsch PetscCall(TSSetTime(ts, ctx->tfaulton));
10989566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, ctx->tfaultoff));
10999566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, X));
1100c4762a1bSJed Brown
1101c4762a1bSJed Brown /* Remove the fault */
11029371c9d4SSatish Balay row_loc = 2 * ctx->faultbus;
11039371c9d4SSatish Balay col_loc = 2 * ctx->faultbus + 1;
1104c4762a1bSJed Brown val = -1 / ctx->Rfault;
11059566063dSJacob Faibussowitsch PetscCall(MatSetValues(ctx->Ybus, 1, &row_loc, 1, &col_loc, &val, ADD_VALUES));
11069371c9d4SSatish Balay row_loc = 2 * ctx->faultbus + 1;
11079371c9d4SSatish Balay col_loc = 2 * ctx->faultbus;
1108c4762a1bSJed Brown val = -1 / ctx->Rfault;
11099566063dSJacob Faibussowitsch PetscCall(MatSetValues(ctx->Ybus, 1, &row_loc, 1, &col_loc, &val, ADD_VALUES));
1110c4762a1bSJed Brown
11119566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(ctx->Ybus, MAT_FINAL_ASSEMBLY));
11129566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(ctx->Ybus, MAT_FINAL_ASSEMBLY));
1113c4762a1bSJed Brown
11149566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(J));
1115c4762a1bSJed Brown
1116c4762a1bSJed Brown ctx->alg_flg = PETSC_TRUE;
1117c4762a1bSJed Brown
1118c4762a1bSJed Brown /* Solve the algebraic equations */
11199566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes_alg, NULL, X));
1120c4762a1bSJed Brown
1121c4762a1bSJed Brown ctx->stepnum++;
1122c4762a1bSJed Brown
1123c4762a1bSJed Brown /* Post-disturbance period */
1124c4762a1bSJed Brown ctx->alg_flg = PETSC_TRUE;
11259566063dSJacob Faibussowitsch PetscCall(TSSetTime(ts, ctx->tfaultoff));
11269566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, ctx->tmax));
11279566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, X));
1128c4762a1bSJed Brown
11299566063dSJacob Faibussowitsch PetscCall(VecGetArray(ctx->vec_q, &x_ptr));
1130c4762a1bSJed Brown *f = x_ptr[0];
11319566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(ctx->vec_q, &x_ptr));
1132c4762a1bSJed Brown
11339566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ctx->Ybus));
11349566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx->V0));
11359566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes_alg));
11369566063dSJacob Faibussowitsch PetscCall(VecDestroy(&F_alg));
11379566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J));
11389566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X));
11399566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts));
1140c4762a1bSJed Brown
1141c4762a1bSJed Brown return 0;
1142c4762a1bSJed Brown }
1143c4762a1bSJed Brown
1144c4762a1bSJed Brown /*TEST
1145c4762a1bSJed Brown
1146c4762a1bSJed Brown build:
1147dfd57a17SPierre Jolivet requires: double !complex !defined(USE_64BIT_INDICES)
1148c4762a1bSJed Brown
1149c4762a1bSJed Brown test:
1150c4762a1bSJed Brown args: -viewer_binary_skip_info -tao_monitor -tao_gttol .2
1151c4762a1bSJed Brown localrunfiles: petscoptions X.bin Ybus.bin
1152c4762a1bSJed Brown
1153c4762a1bSJed Brown TEST*/
1154