Actual source code: ex31.c
slepc-main 2024-11-09
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
11: static char help[] = "Power grid small signal stability analysis of WECC 9 bus system.\n\
12: This example is based on the 9-bus (node) example given in the book Power\n\
13: Systems Dynamics and Stability (Chapter 8) by P. Sauer and M. A. Pai.\n\
14: The power grid in this example consists of 9 buses (nodes), 3 generators,\n\
15: 3 loads, and 9 transmission lines. The network equations are written\n\
16: in current balance form using rectangular coordinates. It uses the SLEPc\n\
17: package to calculate the eigenvalues for small signal stability analysis\n\n";
19: /*
20: This example is based on PETSc's ex9bus example (under TS).
22: The equations for the stability analysis are described by the DAE
24: \dot{x} = f(x,y,t)
25: 0 = g(x,y,t)
27: where the generators are described by differential equations, while the algebraic
28: constraints define the network equations.
30: The generators are modeled with a 4th order differential equation describing the electrical
31: and mechanical dynamics. Each generator also has an exciter system modeled by 3rd order
32: diff. eqns. describing the exciter, voltage regulator, and the feedback stabilizer
33: mechanism.
35: The network equations are described by nodal current balance equations.
36: I(x,y) - Y*V = 0
38: where:
39: I(x,y) is the current injected from generators and loads.
40: Y is the admittance matrix, and
41: V is the voltage vector
43: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
45: The linearized equations for the eigenvalue analysis are
47: \dot{\delta{x}} = f_x*\delta{x} + f_y*\delta{y}
48: 0 = g_x*\delta{x} + g_y*\delta{y}
50: This gives the linearized sensitivity matrix
51: A = | f_x f_y |
52: | g_x g_y |
54: We are interested in the eigenvalues of the Schur complement of A
55: \hat{A} = f_x - g_x*inv(g_y)*f_y
57: Example contributed by: Shrirang Abhyankar
58: */
60: #include <petscdm.h>
61: #include <petscdmda.h>
62: #include <petscdmcomposite.h>
63: #include <slepceps.h>
65: #define freq 60
66: #define w_s (2*PETSC_PI*freq)
68: /* Sizes and indices */
69: const PetscInt nbus = 9; /* Number of network buses */
70: const PetscInt ngen = 3; /* Number of generators */
71: const PetscInt nload = 3; /* Number of loads */
72: const PetscInt gbus[3] = {0,1,2}; /* Buses at which generators are incident */
73: const PetscInt lbus[3] = {4,5,7}; /* Buses at which loads are incident */
75: /* Generator real and reactive powers (found via loadflow) */
76: const PetscScalar PG[3] = {0.716786142395021,1.630000000000000,0.850000000000000};
77: const PetscScalar QG[3] = {0.270702180178785,0.066120127797275,-0.108402221791588};
78: /* Generator constants */
79: const PetscScalar H[3] = {23.64,6.4,3.01}; /* Inertia constant */
80: const PetscScalar Rs[3] = {0.0,0.0,0.0}; /* Stator Resistance */
81: const PetscScalar Xd[3] = {0.146,0.8958,1.3125}; /* d-axis reactance */
82: const PetscScalar Xdp[3] = {0.0608,0.1198,0.1813}; /* d-axis transient reactance */
83: const PetscScalar Xq[3] = {0.0969,0.8645,1.2578}; /* q-axis reactance Xq(1) set to 0.4360, value given in text 0.0969 */
84: const PetscScalar Xqp[3] = {0.0969,0.1969,0.25}; /* q-axis transient reactance */
85: const PetscScalar Td0p[3] = {8.96,6.0,5.89}; /* d-axis open circuit time constant */
86: const PetscScalar Tq0p[3] = {0.31,0.535,0.6}; /* q-axis open circuit time constant */
87: PetscScalar M[3]; /* M = 2*H/w_s */
88: PetscScalar D[3]; /* D = 0.1*M */
90: PetscScalar TM[3]; /* Mechanical Torque */
91: /* Exciter system constants */
92: const PetscScalar KA[3] = {20.0,20.0,20.0}; /* Voltage regulartor gain constant */
93: const PetscScalar TA[3] = {0.2,0.2,0.2}; /* Voltage regulator time constant */
94: const PetscScalar KE[3] = {1.0,1.0,1.0}; /* Exciter gain constant */
95: const PetscScalar TE[3] = {0.314,0.314,0.314}; /* Exciter time constant */
96: const PetscScalar KF[3] = {0.063,0.063,0.063}; /* Feedback stabilizer gain constant */
97: const PetscScalar TF[3] = {0.35,0.35,0.35}; /* Feedback stabilizer time constant */
98: const PetscScalar k1[3] = {0.0039,0.0039,0.0039};
99: const PetscScalar k2[3] = {1.555,1.555,1.555}; /* k1 and k2 for calculating the saturation function SE = k1*exp(k2*Efd) */
101: PetscScalar Vref[3];
102: /* Load constants
103: We use a composite load model that describes the load and reactive powers at each time instant as follows
104: P(t) = \sum\limits_{i=0}^ld_nsegsp \ld_alphap_i*P_D0(\frac{V_m(t)}{V_m0})^\ld_betap_i
105: Q(t) = \sum\limits_{i=0}^ld_nsegsq \ld_alphaq_i*Q_D0(\frac{V_m(t)}{V_m0})^\ld_betaq_i
106: where
107: ld_nsegsp,ld_nsegsq - Number of individual load models for real and reactive power loads
108: ld_alphap,ld_alphap - Percentage contribution (weights) or loads
109: P_D0 - Real power load
110: Q_D0 - Reactive power load
111: V_m(t) - Voltage magnitude at time t
112: V_m0 - Voltage magnitude at t = 0
113: ld_betap, ld_betaq - exponents describing the load model for real and reactive part
115: Note: All loads have the same characteristic currently.
116: */
117: const PetscScalar PD0[3] = {1.25,0.9,1.0};
118: const PetscScalar QD0[3] = {0.5,0.3,0.35};
119: const PetscInt ld_nsegsp[3] = {3,3,3};
120: const PetscScalar ld_alphap[3] = {0.0,0.0,1.0};
121: const PetscScalar ld_betap[3] = {2.0,1.0,0.0};
122: const PetscInt ld_nsegsq[3] = {3,3,3};
123: const PetscScalar ld_alphaq[3] = {0.0,0.0,1.0};
124: const PetscScalar ld_betaq[3] = {2.0,1.0,0.0};
126: typedef struct {
127: DM dmgen, dmnet; /* DMs to manage generator and network subsystem */
128: DM dmpgrid; /* Composite DM to manage the entire power grid */
129: Mat Ybus; /* Network admittance matrix */
130: Vec V0; /* Initial voltage vector (Power flow solution) */
131: PetscInt neqs_gen,neqs_net,neqs_pgrid;
132: IS is_diff; /* indices for differential equations */
133: IS is_alg; /* indices for algebraic equations */
134: } Userctx;
136: /* Converts from machine frame (dq) to network (phase a real,imag) reference frame */
137: PetscErrorCode dq2ri(PetscScalar Fd,PetscScalar Fq,PetscScalar delta,PetscScalar *Fr,PetscScalar *Fi)
138: {
139: PetscFunctionBegin;
140: *Fr = Fd*PetscSinScalar(delta) + Fq*PetscCosScalar(delta);
141: *Fi = -Fd*PetscCosScalar(delta) + Fq*PetscSinScalar(delta);
142: PetscFunctionReturn(PETSC_SUCCESS);
143: }
145: /* Converts from network frame ([phase a real,imag) to machine (dq) reference frame */
146: PetscErrorCode ri2dq(PetscScalar Fr,PetscScalar Fi,PetscScalar delta,PetscScalar *Fd,PetscScalar *Fq)
147: {
148: PetscFunctionBegin;
149: *Fd = Fr*PetscSinScalar(delta) - Fi*PetscCosScalar(delta);
150: *Fq = Fr*PetscCosScalar(delta) + Fi*PetscSinScalar(delta);
151: PetscFunctionReturn(PETSC_SUCCESS);
152: }
154: PetscErrorCode SetInitialGuess(Vec X,Userctx *user)
155: {
156: Vec Xgen,Xnet;
157: PetscScalar *xgen,*xnet;
158: PetscInt i,idx=0;
159: PetscScalar Vr,Vi,IGr,IGi,Vm,Vm2;
160: PetscScalar Eqp,Edp,delta;
161: PetscScalar Efd,RF,VR; /* Exciter variables */
162: PetscScalar Id,Iq; /* Generator dq axis currents */
163: PetscScalar theta,Vd,Vq,SE;
165: PetscFunctionBegin;
166: M[0] = 2*H[0]/w_s; M[1] = 2*H[1]/w_s; M[2] = 2*H[2]/w_s;
167: /* D[0] = 0.1*M[0]; D[1] = 0.1*M[1]; D[2] = 0.1*M[2];
168: */
169: D[0] = D[1] = D[2] = 0.0;
170: PetscCall(DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet));
172: /* Network subsystem initialization */
173: PetscCall(VecCopy(user->V0,Xnet));
175: /* Generator subsystem initialization */
176: PetscCall(VecGetArray(Xgen,&xgen));
177: PetscCall(VecGetArray(Xnet,&xnet));
179: for (i=0; i < ngen; i++) {
180: Vr = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */
181: Vi = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */
182: Vm = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm;
183: IGr = (Vr*PG[i] + Vi*QG[i])/Vm2;
184: IGi = (Vi*PG[i] - Vr*QG[i])/Vm2;
186: delta = atan2(Vi+Xq[i]*IGr,Vr-Xq[i]*IGi); /* Machine angle */
188: theta = PETSC_PI/2.0 - delta;
190: Id = IGr*PetscCosScalar(theta) - IGi*PetscSinScalar(theta); /* d-axis stator current */
191: Iq = IGr*PetscSinScalar(theta) + IGi*PetscCosScalar(theta); /* q-axis stator current */
193: Vd = Vr*PetscCosScalar(theta) - Vi*PetscSinScalar(theta);
194: Vq = Vr*PetscSinScalar(theta) + Vi*PetscCosScalar(theta);
196: Edp = Vd + Rs[i]*Id - Xqp[i]*Iq; /* d-axis transient EMF */
197: Eqp = Vq + Rs[i]*Iq + Xdp[i]*Id; /* q-axis transient EMF */
199: TM[i] = PG[i];
201: /* The generator variables are ordered as [Eqp,Edp,delta,w,Id,Iq] */
202: xgen[idx] = Eqp;
203: xgen[idx+1] = Edp;
204: xgen[idx+2] = delta;
205: xgen[idx+3] = w_s;
207: idx = idx + 4;
209: xgen[idx] = Id;
210: xgen[idx+1] = Iq;
212: idx = idx + 2;
214: /* Exciter */
215: Efd = Eqp + (Xd[i] - Xdp[i])*Id;
216: SE = k1[i]*PetscExpScalar(k2[i]*Efd);
217: VR = KE[i]*Efd + SE;
218: RF = KF[i]*Efd/TF[i];
220: xgen[idx] = Efd;
221: xgen[idx+1] = RF;
222: xgen[idx+2] = VR;
224: Vref[i] = Vm + (VR/KA[i]);
226: idx = idx + 3;
227: }
229: PetscCall(VecRestoreArray(Xgen,&xgen));
230: PetscCall(VecRestoreArray(Xnet,&xnet));
232: /* PetscCall(VecView(Xgen,0)); */
233: PetscCall(DMCompositeGather(user->dmpgrid,INSERT_VALUES,X,Xgen,Xnet));
234: PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet));
235: PetscFunctionReturn(PETSC_SUCCESS);
236: }
238: PetscErrorCode PreallocateJacobian(Mat J,Userctx *user)
239: {
240: PetscInt *d_nnz;
241: PetscInt i,idx=0,start=0;
242: PetscInt ncols;
244: PetscFunctionBegin;
245: PetscCall(PetscMalloc1(user->neqs_pgrid,&d_nnz));
246: for (i=0; i<user->neqs_pgrid; i++) d_nnz[i] = 0;
247: /* Generator subsystem */
248: for (i=0; i < ngen; i++) {
250: d_nnz[idx] += 3;
251: d_nnz[idx+1] += 2;
252: d_nnz[idx+2] += 2;
253: d_nnz[idx+3] += 5;
254: d_nnz[idx+4] += 6;
255: d_nnz[idx+5] += 6;
257: d_nnz[user->neqs_gen+2*gbus[i]] += 3;
258: d_nnz[user->neqs_gen+2*gbus[i]+1] += 3;
260: d_nnz[idx+6] += 2;
261: d_nnz[idx+7] += 2;
262: d_nnz[idx+8] += 5;
264: idx = idx + 9;
265: }
267: start = user->neqs_gen;
269: for (i=0; i < nbus; i++) {
270: PetscCall(MatGetRow(user->Ybus,2*i,&ncols,NULL,NULL));
271: d_nnz[start+2*i] += ncols;
272: d_nnz[start+2*i+1] += ncols;
273: PetscCall(MatRestoreRow(user->Ybus,2*i,&ncols,NULL,NULL));
274: }
276: PetscCall(MatSeqAIJSetPreallocation(J,0,d_nnz));
278: PetscCall(PetscFree(d_nnz));
279: PetscFunctionReturn(PETSC_SUCCESS);
280: }
282: /*
283: J = [-df_dx, -df_dy
284: dg_dx, dg_dy]
285: */
286: PetscErrorCode ResidualJacobian(Vec X,Mat J,void *ctx)
287: {
288: Userctx *user=(Userctx*)ctx;
289: Vec Xgen,Xnet;
290: PetscScalar *xgen,*xnet;
291: PetscInt i,idx=0;
292: PetscScalar Vr,Vi,Vm,Vm2;
293: PetscScalar Eqp,Edp,delta; /* Generator variables */
294: PetscScalar Efd;
295: PetscScalar Id,Iq; /* Generator dq axis currents */
296: PetscScalar Vd,Vq;
297: PetscScalar val[10];
298: PetscInt row[2],col[10];
299: PetscInt net_start=user->neqs_gen;
300: PetscScalar Zdq_inv[4],det;
301: PetscScalar dVd_dVr,dVd_dVi,dVq_dVr,dVq_dVi,dVd_ddelta,dVq_ddelta;
302: PetscScalar dIGr_ddelta,dIGi_ddelta,dIGr_dId,dIGr_dIq,dIGi_dId,dIGi_dIq;
303: PetscScalar dSE_dEfd;
304: PetscScalar dVm_dVd,dVm_dVq,dVm_dVr,dVm_dVi;
305: PetscInt ncols;
306: const PetscInt *cols;
307: const PetscScalar *yvals;
308: PetscInt k;
309: PetscScalar PD,QD,Vm0,*v0,Vm4;
310: PetscScalar dPD_dVr,dPD_dVi,dQD_dVr,dQD_dVi;
311: PetscScalar dIDr_dVr,dIDr_dVi,dIDi_dVr,dIDi_dVi;
313: PetscFunctionBegin;
314: PetscCall(MatZeroEntries(J));
315: PetscCall(DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet));
316: PetscCall(DMCompositeScatter(user->dmpgrid,X,Xgen,Xnet));
318: PetscCall(VecGetArray(Xgen,&xgen));
319: PetscCall(VecGetArray(Xnet,&xnet));
321: /* Generator subsystem */
322: for (i=0; i < ngen; i++) {
323: Eqp = xgen[idx];
324: Edp = xgen[idx+1];
325: delta = xgen[idx+2];
326: Id = xgen[idx+4];
327: Iq = xgen[idx+5];
328: Efd = xgen[idx+6];
330: /* fgen[idx] = (Eqp + (Xd[i] - Xdp[i])*Id - Efd)/Td0p[i]; */
331: row[0] = idx;
332: col[0] = idx; col[1] = idx+4; col[2] = idx+6;
333: val[0] = 1/ Td0p[i]; val[1] = (Xd[i] - Xdp[i])/ Td0p[i]; val[2] = -1/Td0p[i];
335: PetscCall(MatSetValues(J,1,row,3,col,val,INSERT_VALUES));
337: /* fgen[idx+1] = (Edp - (Xq[i] - Xqp[i])*Iq)/Tq0p[i]; */
338: row[0] = idx + 1;
339: col[0] = idx + 1; col[1] = idx+5;
340: val[0] = 1/Tq0p[i]; val[1] = -(Xq[i] - Xqp[i])/Tq0p[i];
341: PetscCall(MatSetValues(J,1,row,2,col,val,INSERT_VALUES));
343: /* fgen[idx+2] = - w + w_s; */
344: row[0] = idx + 2;
345: col[0] = idx + 2; col[1] = idx + 3;
346: val[0] = 0; val[1] = -1;
347: PetscCall(MatSetValues(J,1,row,2,col,val,INSERT_VALUES));
349: /* fgen[idx+3] = (-TM[i] + Edp*Id + Eqp*Iq + (Xqp[i] - Xdp[i])*Id*Iq + D[i]*(w - w_s))/M[i]; */
350: row[0] = idx + 3;
351: col[0] = idx; col[1] = idx + 1; col[2] = idx + 3; col[3] = idx + 4; col[4] = idx + 5;
352: val[0] = Iq/M[i]; val[1] = Id/M[i]; val[2] = D[i]/M[i]; val[3] = (Edp + (Xqp[i]-Xdp[i])*Iq)/M[i]; val[4] = (Eqp + (Xqp[i] - Xdp[i])*Id)/M[i];
353: PetscCall(MatSetValues(J,1,row,5,col,val,INSERT_VALUES));
355: Vr = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */
356: Vi = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */
357: PetscCall(ri2dq(Vr,Vi,delta,&Vd,&Vq));
359: det = Rs[i]*Rs[i] + Xdp[i]*Xqp[i];
361: Zdq_inv[0] = Rs[i]/det;
362: Zdq_inv[1] = Xqp[i]/det;
363: Zdq_inv[2] = -Xdp[i]/det;
364: Zdq_inv[3] = Rs[i]/det;
366: dVd_dVr = PetscSinScalar(delta); dVd_dVi = -PetscCosScalar(delta);
367: dVq_dVr = PetscCosScalar(delta); dVq_dVi = PetscSinScalar(delta);
368: dVd_ddelta = Vr*PetscCosScalar(delta) + Vi*PetscSinScalar(delta);
369: dVq_ddelta = -Vr*PetscSinScalar(delta) + Vi*PetscCosScalar(delta);
371: /* fgen[idx+4] = Zdq_inv[0]*(-Edp + Vd) + Zdq_inv[1]*(-Eqp + Vq) + Id; */
372: row[0] = idx+4;
373: col[0] = idx; col[1] = idx+1; col[2] = idx + 2;
374: val[0] = -Zdq_inv[1]; val[1] = -Zdq_inv[0]; val[2] = Zdq_inv[0]*dVd_ddelta + Zdq_inv[1]*dVq_ddelta;
375: col[3] = idx + 4; col[4] = net_start+2*gbus[i]; col[5] = net_start + 2*gbus[i]+1;
376: val[3] = 1; val[4] = Zdq_inv[0]*dVd_dVr + Zdq_inv[1]*dVq_dVr; val[5] = Zdq_inv[0]*dVd_dVi + Zdq_inv[1]*dVq_dVi;
377: PetscCall(MatSetValues(J,1,row,6,col,val,INSERT_VALUES));
379: /* fgen[idx+5] = Zdq_inv[2]*(-Edp + Vd) + Zdq_inv[3]*(-Eqp + Vq) + Iq; */
380: row[0] = idx+5;
381: col[0] = idx; col[1] = idx+1; col[2] = idx + 2;
382: val[0] = -Zdq_inv[3]; val[1] = -Zdq_inv[2]; val[2] = Zdq_inv[2]*dVd_ddelta + Zdq_inv[3]*dVq_ddelta;
383: col[3] = idx + 5; col[4] = net_start+2*gbus[i]; col[5] = net_start + 2*gbus[i]+1;
384: val[3] = 1; val[4] = Zdq_inv[2]*dVd_dVr + Zdq_inv[3]*dVq_dVr; val[5] = Zdq_inv[2]*dVd_dVi + Zdq_inv[3]*dVq_dVi;
385: PetscCall(MatSetValues(J,1,row,6,col,val,INSERT_VALUES));
387: dIGr_ddelta = Id*PetscCosScalar(delta) - Iq*PetscSinScalar(delta);
388: dIGi_ddelta = Id*PetscSinScalar(delta) + Iq*PetscCosScalar(delta);
389: dIGr_dId = PetscSinScalar(delta); dIGr_dIq = PetscCosScalar(delta);
390: dIGi_dId = -PetscCosScalar(delta); dIGi_dIq = PetscSinScalar(delta);
392: /* fnet[2*gbus[i]] -= IGi; */
393: row[0] = net_start + 2*gbus[i];
394: col[0] = idx+2; col[1] = idx + 4; col[2] = idx + 5;
395: val[0] = -dIGi_ddelta; val[1] = -dIGi_dId; val[2] = -dIGi_dIq;
396: PetscCall(MatSetValues(J,1,row,3,col,val,INSERT_VALUES));
398: /* fnet[2*gbus[i]+1] -= IGr; */
399: row[0] = net_start + 2*gbus[i]+1;
400: col[0] = idx+2; col[1] = idx + 4; col[2] = idx + 5;
401: val[0] = -dIGr_ddelta; val[1] = -dIGr_dId; val[2] = -dIGr_dIq;
402: PetscCall(MatSetValues(J,1,row,3,col,val,INSERT_VALUES));
404: Vm = PetscSqrtScalar(Vd*Vd + Vq*Vq); Vm2 = Vm*Vm;
406: /* fgen[idx+6] = (KE[i]*Efd + SE - VR)/TE[i]; */
407: /* SE = k1[i]*PetscExpScalar(k2[i]*Efd); */
409: dSE_dEfd = k1[i]*k2[i]*PetscExpScalar(k2[i]*Efd);
411: row[0] = idx + 6;
412: col[0] = idx + 6; col[1] = idx + 8;
413: val[0] = (KE[i] + dSE_dEfd)/TE[i]; val[1] = -1/TE[i];
414: PetscCall(MatSetValues(J,1,row,2,col,val,INSERT_VALUES));
416: /* Exciter differential equations */
418: /* fgen[idx+7] = (RF - KF[i]*Efd/TF[i])/TF[i]; */
419: row[0] = idx + 7;
420: col[0] = idx + 6; col[1] = idx + 7;
421: val[0] = (-KF[i]/TF[i])/TF[i]; val[1] = 1/TF[i];
422: PetscCall(MatSetValues(J,1,row,2,col,val,INSERT_VALUES));
424: /* fgen[idx+8] = (VR - KA[i]*RF + KA[i]*KF[i]*Efd/TF[i] - KA[i]*(Vref[i] - Vm))/TA[i]; */
425: /* Vm = (Vd^2 + Vq^2)^0.5; */
427: dVm_dVd = Vd/Vm; dVm_dVq = Vq/Vm;
428: dVm_dVr = dVm_dVd*dVd_dVr + dVm_dVq*dVq_dVr;
429: dVm_dVi = dVm_dVd*dVd_dVi + dVm_dVq*dVq_dVi;
430: row[0] = idx + 8;
431: col[0] = idx + 6; col[1] = idx + 7; col[2] = idx + 8;
432: val[0] = (KA[i]*KF[i]/TF[i])/TA[i]; val[1] = -KA[i]/TA[i]; val[2] = 1/TA[i];
433: col[3] = net_start + 2*gbus[i]; col[4] = net_start + 2*gbus[i]+1;
434: val[3] = KA[i]*dVm_dVr/TA[i]; val[4] = KA[i]*dVm_dVi/TA[i];
435: PetscCall(MatSetValues(J,1,row,5,col,val,INSERT_VALUES));
436: idx = idx + 9;
437: }
439: for (i=0; i<nbus; i++) {
440: PetscCall(MatGetRow(user->Ybus,2*i,&ncols,&cols,&yvals));
441: row[0] = net_start + 2*i;
442: for (k=0; k<ncols; k++) {
443: col[k] = net_start + cols[k];
444: val[k] = yvals[k];
445: }
446: PetscCall(MatSetValues(J,1,row,ncols,col,val,INSERT_VALUES));
447: PetscCall(MatRestoreRow(user->Ybus,2*i,&ncols,&cols,&yvals));
449: PetscCall(MatGetRow(user->Ybus,2*i+1,&ncols,&cols,&yvals));
450: row[0] = net_start + 2*i+1;
451: for (k=0; k<ncols; k++) {
452: col[k] = net_start + cols[k];
453: val[k] = yvals[k];
454: }
455: PetscCall(MatSetValues(J,1,row,ncols,col,val,INSERT_VALUES));
456: PetscCall(MatRestoreRow(user->Ybus,2*i+1,&ncols,&cols,&yvals));
457: }
459: PetscCall(MatAssemblyBegin(J,MAT_FLUSH_ASSEMBLY));
460: PetscCall(MatAssemblyEnd(J,MAT_FLUSH_ASSEMBLY));
462: PetscCall(VecGetArray(user->V0,&v0));
463: for (i=0; i < nload; i++) {
464: Vr = xnet[2*lbus[i]]; /* Real part of load bus voltage */
465: Vi = xnet[2*lbus[i]+1]; /* Imaginary part of the load bus voltage */
466: Vm = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm; Vm4 = Vm2*Vm2;
467: Vm0 = PetscSqrtScalar(v0[2*lbus[i]]*v0[2*lbus[i]] + v0[2*lbus[i]+1]*v0[2*lbus[i]+1]);
468: PD = QD = 0.0;
469: dPD_dVr = dPD_dVi = dQD_dVr = dQD_dVi = 0.0;
470: for (k=0; k < ld_nsegsp[i]; k++) {
471: PD += ld_alphap[k]*PD0[i]*PetscPowScalar((Vm/Vm0),ld_betap[k]);
472: dPD_dVr += ld_alphap[k]*ld_betap[k]*PD0[i]*PetscPowScalar((1/Vm0),ld_betap[k])*Vr*PetscPowScalar(Vm,(ld_betap[k]-2));
473: dPD_dVi += ld_alphap[k]*ld_betap[k]*PD0[i]*PetscPowScalar((1/Vm0),ld_betap[k])*Vi*PetscPowScalar(Vm,(ld_betap[k]-2));
474: }
475: for (k=0; k < ld_nsegsq[i]; k++) {
476: QD += ld_alphaq[k]*QD0[i]*PetscPowScalar((Vm/Vm0),ld_betaq[k]);
477: dQD_dVr += ld_alphaq[k]*ld_betaq[k]*QD0[i]*PetscPowScalar((1/Vm0),ld_betaq[k])*Vr*PetscPowScalar(Vm,(ld_betaq[k]-2));
478: dQD_dVi += ld_alphaq[k]*ld_betaq[k]*QD0[i]*PetscPowScalar((1/Vm0),ld_betaq[k])*Vi*PetscPowScalar(Vm,(ld_betaq[k]-2));
479: }
481: /* IDr = (PD*Vr + QD*Vi)/Vm2; */
482: /* IDi = (-QD*Vr + PD*Vi)/Vm2; */
484: dIDr_dVr = (dPD_dVr*Vr + dQD_dVr*Vi + PD)/Vm2 - ((PD*Vr + QD*Vi)*2*Vr)/Vm4;
485: dIDr_dVi = (dPD_dVi*Vr + dQD_dVi*Vi + QD)/Vm2 - ((PD*Vr + QD*Vi)*2*Vi)/Vm4;
487: dIDi_dVr = (-dQD_dVr*Vr + dPD_dVr*Vi - QD)/Vm2 - ((-QD*Vr + PD*Vi)*2*Vr)/Vm4;
488: dIDi_dVi = (-dQD_dVi*Vr + dPD_dVi*Vi + PD)/Vm2 - ((-QD*Vr + PD*Vi)*2*Vi)/Vm4;
490: /* fnet[2*lbus[i]] += IDi; */
491: row[0] = net_start + 2*lbus[i];
492: col[0] = net_start + 2*lbus[i]; col[1] = net_start + 2*lbus[i]+1;
493: val[0] = dIDi_dVr; val[1] = dIDi_dVi;
494: PetscCall(MatSetValues(J,1,row,2,col,val,ADD_VALUES));
495: /* fnet[2*lbus[i]+1] += IDr; */
496: row[0] = net_start + 2*lbus[i]+1;
497: col[0] = net_start + 2*lbus[i]; col[1] = net_start + 2*lbus[i]+1;
498: val[0] = dIDr_dVr; val[1] = dIDr_dVi;
499: PetscCall(MatSetValues(J,1,row,2,col,val,ADD_VALUES));
500: }
501: PetscCall(VecRestoreArray(user->V0,&v0));
503: PetscCall(VecRestoreArray(Xgen,&xgen));
504: PetscCall(VecRestoreArray(Xnet,&xnet));
506: PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet));
508: PetscCall(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY));
509: PetscCall(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY));
510: PetscFunctionReturn(PETSC_SUCCESS);
511: }
513: int main(int argc,char **argv)
514: {
515: EPS eps;
516: EPSType type;
517: PetscMPIInt size;
518: Userctx user;
519: PetscViewer Xview,Ybusview;
520: Vec X,Xr,Xi;
521: Mat J,Jred=NULL;
522: IS is0,is1;
523: PetscInt i,*idx2,its,nev,nconv;
524: PetscReal error,re,im;
525: PetscScalar kr,ki;
526: PetscBool terse;
528: PetscFunctionBeginUser;
529: PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
530: PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
531: PetscCheck(size==1,PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs");
532: /* show detailed info unless -terse option is given by user */
533: PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
535: user.neqs_gen = 9*ngen; /* # eqs. for generator subsystem */
536: user.neqs_net = 2*nbus; /* # eqs. for network subsystem */
537: user.neqs_pgrid = user.neqs_gen + user.neqs_net;
538: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nStability analysis in a network with %" PetscInt_FMT " buses and %" PetscInt_FMT " generators\n\n",nbus,ngen));
540: /* Create indices for differential and algebraic equations */
541: PetscCall(PetscMalloc1(7*ngen,&idx2));
542: for (i=0; i<ngen; i++) {
543: idx2[7*i] = 9*i; idx2[7*i+1] = 9*i+1; idx2[7*i+2] = 9*i+2; idx2[7*i+3] = 9*i+3;
544: idx2[7*i+4] = 9*i+6; idx2[7*i+5] = 9*i+7; idx2[7*i+6] = 9*i+8;
545: }
546: PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,7*ngen,idx2,PETSC_COPY_VALUES,&user.is_diff));
547: PetscCall(ISComplement(user.is_diff,0,user.neqs_pgrid,&user.is_alg));
548: PetscCall(PetscFree(idx2));
550: /* Read initial voltage vector and Ybus */
551: PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"X.bin",FILE_MODE_READ,&Xview));
552: PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"Ybus.bin",FILE_MODE_READ,&Ybusview));
554: PetscCall(VecCreate(PETSC_COMM_WORLD,&user.V0));
555: PetscCall(VecSetSizes(user.V0,PETSC_DECIDE,user.neqs_net));
556: PetscCall(VecLoad(user.V0,Xview));
558: PetscCall(MatCreate(PETSC_COMM_WORLD,&user.Ybus));
559: PetscCall(MatSetSizes(user.Ybus,PETSC_DECIDE,PETSC_DECIDE,user.neqs_net,user.neqs_net));
560: PetscCall(MatSetType(user.Ybus,MATBAIJ));
561: /* PetscCall(MatSetBlockSize(user.Ybus,2)); */
562: PetscCall(MatLoad(user.Ybus,Ybusview));
564: PetscCall(PetscViewerDestroy(&Xview));
565: PetscCall(PetscViewerDestroy(&Ybusview));
567: /* Create DMs for generator and network subsystems */
568: PetscCall(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,user.neqs_gen,1,1,NULL,&user.dmgen));
569: PetscCall(DMSetOptionsPrefix(user.dmgen,"dmgen_"));
570: PetscCall(DMSetFromOptions(user.dmgen));
571: PetscCall(DMSetUp(user.dmgen));
572: PetscCall(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,user.neqs_net,1,1,NULL,&user.dmnet));
573: PetscCall(DMSetOptionsPrefix(user.dmnet,"dmnet_"));
574: PetscCall(DMSetFromOptions(user.dmnet));
575: PetscCall(DMSetUp(user.dmnet));
577: /* Create a composite DM packer and add the two DMs */
578: PetscCall(DMCompositeCreate(PETSC_COMM_WORLD,&user.dmpgrid));
579: PetscCall(DMSetOptionsPrefix(user.dmpgrid,"pgrid_"));
580: PetscCall(DMCompositeAddDM(user.dmpgrid,user.dmgen));
581: PetscCall(DMCompositeAddDM(user.dmpgrid,user.dmnet));
583: PetscCall(DMCreateGlobalVector(user.dmpgrid,&X));
585: PetscCall(MatCreate(PETSC_COMM_WORLD,&J));
586: PetscCall(MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,user.neqs_pgrid,user.neqs_pgrid));
587: PetscCall(MatSetFromOptions(J));
588: PetscCall(PreallocateJacobian(J,&user));
590: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
591: Set initial conditions
592: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
593: PetscCall(SetInitialGuess(X,&user));
595: /* Form Jacobian */
596: PetscCall(ResidualJacobian(X,J,(void*)&user));
597: PetscCall(MatScale(J,-1));
598: is0 = user.is_diff;
599: is1 = user.is_alg;
601: PetscCall(MatGetSchurComplement(J,is1,is1,is0,is0,MAT_IGNORE_MATRIX,NULL,MAT_SCHUR_COMPLEMENT_AINV_DIAG,MAT_INITIAL_MATRIX,&Jred));
603: if (!terse) PetscCall(MatView(Jred,NULL));
605: PetscCall(MatCreateVecs(Jred,NULL,&Xr));
606: PetscCall(MatCreateVecs(Jred,NULL,&Xi));
608: /* Create the eigensolver and set the various options */
609: PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
610: PetscCall(EPSSetOperators(eps,Jred,NULL));
611: PetscCall(EPSSetProblemType(eps,EPS_NHEP));
612: PetscCall(EPSSetFromOptions(eps));
614: /* Solve the eigenvalue problem */
615: PetscCall(EPSSolve(eps));
617: PetscCall(EPSGetIterationNumber(eps,&its));
618: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the eigensolver: %" PetscInt_FMT "\n",its));
619: PetscCall(EPSGetType(eps,&type));
620: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n", type));
621: PetscCall(EPSGetDimensions(eps,&nev,NULL,NULL));
622: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev));
624: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
625: Display solution and clean up
626: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
627: if (terse) PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL));
628: else {
629: /* Get number of converged approximate eigenpairs */
630: PetscCall(EPSGetConverged(eps,&nconv));
631: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of converged eigenpairs: %" PetscInt_FMT "\n\n",nconv));
633: if (nconv>0) {
634: /* Display eigenvalues and relative errors */
635: PetscCall(PetscPrintf(PETSC_COMM_WORLD,
636: " k ||Ax-kx||/||kx||\n"
637: " ----------------- ------------------\n"));
639: for (i=0;i<nconv;i++) {
640: /* Get converged eigenpairs: i-th eigenvalue is stored in kr (real part) and
641: ki (imaginary part) */
642: PetscCall(EPSGetEigenpair(eps,i,&kr,&ki,Xr,Xi));
643: /* Compute the relative error associated to each eigenpair */
644: PetscCall(EPSComputeError(eps,i,EPS_ERROR_RELATIVE,&error));
646: #if defined(PETSC_USE_COMPLEX)
647: re = PetscRealPart(kr);
648: im = PetscImaginaryPart(kr);
649: #else
650: re = kr;
651: im = ki;
652: #endif
653: if (im!=0.0) PetscCall(PetscPrintf(PETSC_COMM_WORLD," %9f%+9fi %12g\n",(double)re,(double)im,(double)error));
654: else PetscCall(PetscPrintf(PETSC_COMM_WORLD," %12f %12g\n",(double)re,(double)error));
655: }
656: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n"));
657: }
658: }
660: /* Free work space */
661: PetscCall(EPSDestroy(&eps));
662: PetscCall(MatDestroy(&J));
663: PetscCall(MatDestroy(&Jred));
664: PetscCall(MatDestroy(&user.Ybus));
665: PetscCall(VecDestroy(&X));
666: PetscCall(VecDestroy(&Xr));
667: PetscCall(VecDestroy(&Xi));
668: PetscCall(VecDestroy(&user.V0));
669: PetscCall(DMDestroy(&user.dmgen));
670: PetscCall(DMDestroy(&user.dmnet));
671: PetscCall(DMDestroy(&user.dmpgrid));
672: PetscCall(ISDestroy(&user.is_diff));
673: PetscCall(ISDestroy(&user.is_alg));
674: PetscCall(SlepcFinalize());
675: return 0;
676: }
678: /*TEST
680: build:
681: requires: !complex
683: test:
684: suffix: 1
685: args: -terse
686: requires: double !complex !defined(PETSC_USE_64BIT_INDICES)
687: localrunfiles: X.bin Ybus.bin
689: TEST*/