Actual source code: ex31.c
slepc-3.11.2 2019-07-30
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2019, 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
58: Example contributed by: Shrirang Abhyankar
59: */
61: #include <petscdm.h>
62: #include <petscdmda.h>
63: #include <petscdmcomposite.h>
64: #include <slepceps.h>
66: #define freq 60
67: #define w_s (2*PETSC_PI*freq)
69: /* Sizes and indices */
70: const PetscInt nbus = 9; /* Number of network buses */
71: const PetscInt ngen = 3; /* Number of generators */
72: const PetscInt nload = 3; /* Number of loads */
73: const PetscInt gbus[3] = {0,1,2}; /* Buses at which generators are incident */
74: const PetscInt lbus[3] = {4,5,7}; /* Buses at which loads are incident */
76: /* Generator real and reactive powers (found via loadflow) */
77: const PetscScalar PG[3] = {0.716786142395021,1.630000000000000,0.850000000000000};
78: const PetscScalar QG[3] = {0.270702180178785,0.066120127797275,-0.108402221791588};
79: /* Generator constants */
80: const PetscScalar H[3] = {23.64,6.4,3.01}; /* Inertia constant */
81: const PetscScalar Rs[3] = {0.0,0.0,0.0}; /* Stator Resistance */
82: const PetscScalar Xd[3] = {0.146,0.8958,1.3125}; /* d-axis reactance */
83: const PetscScalar Xdp[3] = {0.0608,0.1198,0.1813}; /* d-axis transient reactance */
84: 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 */
85: const PetscScalar Xqp[3] = {0.0969,0.1969,0.25}; /* q-axis transient reactance */
86: const PetscScalar Td0p[3] = {8.96,6.0,5.89}; /* d-axis open circuit time constant */
87: const PetscScalar Tq0p[3] = {0.31,0.535,0.6}; /* q-axis open circuit time constant */
88: PetscScalar M[3]; /* M = 2*H/w_s */
89: PetscScalar D[3]; /* D = 0.1*M */
91: PetscScalar TM[3]; /* Mechanical Torque */
92: /* Exciter system constants */
93: const PetscScalar KA[3] = {20.0,20.0,20.0}; /* Voltage regulartor gain constant */
94: const PetscScalar TA[3] = {0.2,0.2,0.2}; /* Voltage regulator time constant */
95: const PetscScalar KE[3] = {1.0,1.0,1.0}; /* Exciter gain constant */
96: const PetscScalar TE[3] = {0.314,0.314,0.314}; /* Exciter time constant */
97: const PetscScalar KF[3] = {0.063,0.063,0.063}; /* Feedback stabilizer gain constant */
98: const PetscScalar TF[3] = {0.35,0.35,0.35}; /* Feedback stabilizer time constant */
99: const PetscScalar k1[3] = {0.0039,0.0039,0.0039};
100: const PetscScalar k2[3] = {1.555,1.555,1.555}; /* k1 and k2 for calculating the saturation function SE = k1*exp(k2*Efd) */
102: PetscScalar Vref[3];
103: /* Load constants
104: We use a composite load model that describes the load and reactive powers at each time instant as follows
105: P(t) = \sum\limits_{i=0}^ld_nsegsp \ld_alphap_i*P_D0(\frac{V_m(t)}{V_m0})^\ld_betap_i
106: Q(t) = \sum\limits_{i=0}^ld_nsegsq \ld_alphaq_i*Q_D0(\frac{V_m(t)}{V_m0})^\ld_betaq_i
107: where
108: ld_nsegsp,ld_nsegsq - Number of individual load models for real and reactive power loads
109: ld_alphap,ld_alphap - Percentage contribution (weights) or loads
110: P_D0 - Real power load
111: Q_D0 - Reactive power load
112: V_m(t) - Voltage magnitude at time t
113: V_m0 - Voltage magnitude at t = 0
114: ld_betap, ld_betaq - exponents describing the load model for real and reactive part
116: Note: All loads have the same characteristic currently.
117: */
118: const PetscScalar PD0[3] = {1.25,0.9,1.0};
119: const PetscScalar QD0[3] = {0.5,0.3,0.35};
120: const PetscInt ld_nsegsp[3] = {3,3,3};
121: const PetscScalar ld_alphap[3] = {0.0,0.0,1.0};
122: const PetscScalar ld_betap[3] = {2.0,1.0,0.0};
123: const PetscInt ld_nsegsq[3] = {3,3,3};
124: const PetscScalar ld_alphaq[3] = {0.0,0.0,1.0};
125: const PetscScalar ld_betaq[3] = {2.0,1.0,0.0};
127: typedef struct {
128: DM dmgen, dmnet; /* DMs to manage generator and network subsystem */
129: DM dmpgrid; /* Composite DM to manage the entire power grid */
130: Mat Ybus; /* Network admittance matrix */
131: Vec V0; /* Initial voltage vector (Power flow solution) */
132: PetscInt neqs_gen,neqs_net,neqs_pgrid;
133: IS is_diff; /* indices for differential equations */
134: IS is_alg; /* indices for algebraic equations */
135: } Userctx;
137: /* Converts from machine frame (dq) to network (phase a real,imag) reference frame */
138: PetscErrorCode dq2ri(PetscScalar Fd,PetscScalar Fq,PetscScalar delta,PetscScalar *Fr,PetscScalar *Fi)
139: {
141: *Fr = Fd*PetscSinScalar(delta) + Fq*PetscCosScalar(delta);
142: *Fi = -Fd*PetscCosScalar(delta) + Fq*PetscSinScalar(delta);
143: return(0);
144: }
146: /* Converts from network frame ([phase a real,imag) to machine (dq) reference frame */
147: PetscErrorCode ri2dq(PetscScalar Fr,PetscScalar Fi,PetscScalar delta,PetscScalar *Fd,PetscScalar *Fq)
148: {
150: *Fd = Fr*PetscSinScalar(delta) - Fi*PetscCosScalar(delta);
151: *Fq = Fr*PetscCosScalar(delta) + Fi*PetscSinScalar(delta);
152: return(0);
153: }
155: PetscErrorCode SetInitialGuess(Vec X,Userctx *user)
156: {
158: Vec Xgen,Xnet;
159: PetscScalar *xgen,*xnet;
160: PetscInt i,idx=0;
161: PetscScalar Vr,Vi,IGr,IGi,Vm,Vm2;
162: PetscScalar Eqp,Edp,delta;
163: PetscScalar Efd,RF,VR; /* Exciter variables */
164: PetscScalar Id,Iq; /* Generator dq axis currents */
165: PetscScalar theta,Vd,Vq,SE;
168: M[0] = 2*H[0]/w_s; M[1] = 2*H[1]/w_s; M[2] = 2*H[2]/w_s;
169: /* D[0] = 0.1*M[0]; D[1] = 0.1*M[1]; D[2] = 0.1*M[2];
170: */
171: D[0] = D[1] = D[2] = 0.0;
172: DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
174: /* Network subsystem initialization */
175: VecCopy(user->V0,Xnet);
177: /* Generator subsystem initialization */
178: VecGetArray(Xgen,&xgen);
179: VecGetArray(Xnet,&xnet);
181: for (i=0; i < ngen; i++) {
182: Vr = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */
183: Vi = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */
184: Vm = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm;
185: IGr = (Vr*PG[i] + Vi*QG[i])/Vm2;
186: IGi = (Vi*PG[i] - Vr*QG[i])/Vm2;
188: delta = atan2(Vi+Xq[i]*IGr,Vr-Xq[i]*IGi); /* Machine angle */
190: theta = PETSC_PI/2.0 - delta;
192: Id = IGr*PetscCosScalar(theta) - IGi*PetscSinScalar(theta); /* d-axis stator current */
193: Iq = IGr*PetscSinScalar(theta) + IGi*PetscCosScalar(theta); /* q-axis stator current */
195: Vd = Vr*PetscCosScalar(theta) - Vi*PetscSinScalar(theta);
196: Vq = Vr*PetscSinScalar(theta) + Vi*PetscCosScalar(theta);
198: Edp = Vd + Rs[i]*Id - Xqp[i]*Iq; /* d-axis transient EMF */
199: Eqp = Vq + Rs[i]*Iq + Xdp[i]*Id; /* q-axis transient EMF */
201: TM[i] = PG[i];
203: /* The generator variables are ordered as [Eqp,Edp,delta,w,Id,Iq] */
204: xgen[idx] = Eqp;
205: xgen[idx+1] = Edp;
206: xgen[idx+2] = delta;
207: xgen[idx+3] = w_s;
209: idx = idx + 4;
211: xgen[idx] = Id;
212: xgen[idx+1] = Iq;
214: idx = idx + 2;
216: /* Exciter */
217: Efd = Eqp + (Xd[i] - Xdp[i])*Id;
218: SE = k1[i]*PetscExpScalar(k2[i]*Efd);
219: VR = KE[i]*Efd + SE;
220: RF = KF[i]*Efd/TF[i];
222: xgen[idx] = Efd;
223: xgen[idx+1] = RF;
224: xgen[idx+2] = VR;
226: Vref[i] = Vm + (VR/KA[i]);
228: idx = idx + 3;
229: }
231: VecRestoreArray(Xgen,&xgen);
232: VecRestoreArray(Xnet,&xnet);
234: /* VecView(Xgen,0); */
235: DMCompositeGather(user->dmpgrid,INSERT_VALUES,X,Xgen,Xnet);
236: DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);
237: return(0);
238: }
240: PetscErrorCode PreallocateJacobian(Mat J,Userctx *user)
241: {
243: PetscInt *d_nnz;
244: PetscInt i,idx=0,start=0;
245: PetscInt ncols;
248: PetscMalloc1(user->neqs_pgrid,&d_nnz);
249: for (i=0; i<user->neqs_pgrid; i++) d_nnz[i] = 0;
250: /* Generator subsystem */
251: for (i=0; i < ngen; i++) {
253: d_nnz[idx] += 3;
254: d_nnz[idx+1] += 2;
255: d_nnz[idx+2] += 2;
256: d_nnz[idx+3] += 5;
257: d_nnz[idx+4] += 6;
258: d_nnz[idx+5] += 6;
260: d_nnz[user->neqs_gen+2*gbus[i]] += 3;
261: d_nnz[user->neqs_gen+2*gbus[i]+1] += 3;
263: d_nnz[idx+6] += 2;
264: d_nnz[idx+7] += 2;
265: d_nnz[idx+8] += 5;
267: idx = idx + 9;
268: }
270: start = user->neqs_gen;
272: for (i=0; i < nbus; i++) {
273: MatGetRow(user->Ybus,2*i,&ncols,NULL,NULL);
274: d_nnz[start+2*i] += ncols;
275: d_nnz[start+2*i+1] += ncols;
276: MatRestoreRow(user->Ybus,2*i,&ncols,NULL,NULL);
277: }
279: MatSeqAIJSetPreallocation(J,0,d_nnz);
281: PetscFree(d_nnz);
282: return(0);
283: }
285: /*
286: J = [-df_dx, -df_dy
287: dg_dx, dg_dy]
288: */
289: PetscErrorCode ResidualJacobian(Vec X,Mat J,void *ctx)
290: {
292: Userctx *user=(Userctx*)ctx;
293: Vec Xgen,Xnet;
294: PetscScalar *xgen,*xnet;
295: PetscInt i,idx=0;
296: PetscScalar Vr,Vi,Vm,Vm2;
297: PetscScalar Eqp,Edp,delta; /* Generator variables */
298: PetscScalar Efd;
299: PetscScalar Id,Iq; /* Generator dq axis currents */
300: PetscScalar Vd,Vq;
301: PetscScalar val[10];
302: PetscInt row[2],col[10];
303: PetscInt net_start=user->neqs_gen;
304: PetscScalar Zdq_inv[4],det;
305: PetscScalar dVd_dVr,dVd_dVi,dVq_dVr,dVq_dVi,dVd_ddelta,dVq_ddelta;
306: PetscScalar dIGr_ddelta,dIGi_ddelta,dIGr_dId,dIGr_dIq,dIGi_dId,dIGi_dIq;
307: PetscScalar dSE_dEfd;
308: PetscScalar dVm_dVd,dVm_dVq,dVm_dVr,dVm_dVi;
309: PetscInt ncols;
310: const PetscInt *cols;
311: const PetscScalar *yvals;
312: PetscInt k;
313: PetscScalar PD,QD,Vm0,*v0,Vm4;
314: PetscScalar dPD_dVr,dPD_dVi,dQD_dVr,dQD_dVi;
315: PetscScalar dIDr_dVr,dIDr_dVi,dIDi_dVr,dIDi_dVi;
319: MatZeroEntries(J);
320: DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
321: DMCompositeScatter(user->dmpgrid,X,Xgen,Xnet);
323: VecGetArray(Xgen,&xgen);
324: VecGetArray(Xnet,&xnet);
326: /* Generator subsystem */
327: for (i=0; i < ngen; i++) {
328: Eqp = xgen[idx];
329: Edp = xgen[idx+1];
330: delta = xgen[idx+2];
331: Id = xgen[idx+4];
332: Iq = xgen[idx+5];
333: Efd = xgen[idx+6];
335: /* fgen[idx] = (Eqp + (Xd[i] - Xdp[i])*Id - Efd)/Td0p[i]; */
336: row[0] = idx;
337: col[0] = idx; col[1] = idx+4; col[2] = idx+6;
338: val[0] = 1/ Td0p[i]; val[1] = (Xd[i] - Xdp[i])/ Td0p[i]; val[2] = -1/Td0p[i];
340: MatSetValues(J,1,row,3,col,val,INSERT_VALUES);
342: /* fgen[idx+1] = (Edp - (Xq[i] - Xqp[i])*Iq)/Tq0p[i]; */
343: row[0] = idx + 1;
344: col[0] = idx + 1; col[1] = idx+5;
345: val[0] = 1/Tq0p[i]; val[1] = -(Xq[i] - Xqp[i])/Tq0p[i];
346: MatSetValues(J,1,row,2,col,val,INSERT_VALUES);
348: /* fgen[idx+2] = - w + w_s; */
349: row[0] = idx + 2;
350: col[0] = idx + 2; col[1] = idx + 3;
351: val[0] = 0; val[1] = -1;
352: MatSetValues(J,1,row,2,col,val,INSERT_VALUES);
354: /* fgen[idx+3] = (-TM[i] + Edp*Id + Eqp*Iq + (Xqp[i] - Xdp[i])*Id*Iq + D[i]*(w - w_s))/M[i]; */
355: row[0] = idx + 3;
356: col[0] = idx; col[1] = idx + 1; col[2] = idx + 3; col[3] = idx + 4; col[4] = idx + 5;
357: 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];
358: MatSetValues(J,1,row,5,col,val,INSERT_VALUES);
360: Vr = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */
361: Vi = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */
362: ri2dq(Vr,Vi,delta,&Vd,&Vq);
364: det = Rs[i]*Rs[i] + Xdp[i]*Xqp[i];
366: Zdq_inv[0] = Rs[i]/det;
367: Zdq_inv[1] = Xqp[i]/det;
368: Zdq_inv[2] = -Xdp[i]/det;
369: Zdq_inv[3] = Rs[i]/det;
371: dVd_dVr = PetscSinScalar(delta); dVd_dVi = -PetscCosScalar(delta);
372: dVq_dVr = PetscCosScalar(delta); dVq_dVi = PetscSinScalar(delta);
373: dVd_ddelta = Vr*PetscCosScalar(delta) + Vi*PetscSinScalar(delta);
374: dVq_ddelta = -Vr*PetscSinScalar(delta) + Vi*PetscCosScalar(delta);
376: /* fgen[idx+4] = Zdq_inv[0]*(-Edp + Vd) + Zdq_inv[1]*(-Eqp + Vq) + Id; */
377: row[0] = idx+4;
378: col[0] = idx; col[1] = idx+1; col[2] = idx + 2;
379: val[0] = -Zdq_inv[1]; val[1] = -Zdq_inv[0]; val[2] = Zdq_inv[0]*dVd_ddelta + Zdq_inv[1]*dVq_ddelta;
380: col[3] = idx + 4; col[4] = net_start+2*gbus[i]; col[5] = net_start + 2*gbus[i]+1;
381: 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;
382: MatSetValues(J,1,row,6,col,val,INSERT_VALUES);
384: /* fgen[idx+5] = Zdq_inv[2]*(-Edp + Vd) + Zdq_inv[3]*(-Eqp + Vq) + Iq; */
385: row[0] = idx+5;
386: col[0] = idx; col[1] = idx+1; col[2] = idx + 2;
387: val[0] = -Zdq_inv[3]; val[1] = -Zdq_inv[2]; val[2] = Zdq_inv[2]*dVd_ddelta + Zdq_inv[3]*dVq_ddelta;
388: col[3] = idx + 5; col[4] = net_start+2*gbus[i]; col[5] = net_start + 2*gbus[i]+1;
389: 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;
390: MatSetValues(J,1,row,6,col,val,INSERT_VALUES);
392: dIGr_ddelta = Id*PetscCosScalar(delta) - Iq*PetscSinScalar(delta);
393: dIGi_ddelta = Id*PetscSinScalar(delta) + Iq*PetscCosScalar(delta);
394: dIGr_dId = PetscSinScalar(delta); dIGr_dIq = PetscCosScalar(delta);
395: dIGi_dId = -PetscCosScalar(delta); dIGi_dIq = PetscSinScalar(delta);
397: /* fnet[2*gbus[i]] -= IGi; */
398: row[0] = net_start + 2*gbus[i];
399: col[0] = idx+2; col[1] = idx + 4; col[2] = idx + 5;
400: val[0] = -dIGi_ddelta; val[1] = -dIGi_dId; val[2] = -dIGi_dIq;
401: MatSetValues(J,1,row,3,col,val,INSERT_VALUES);
403: /* fnet[2*gbus[i]+1] -= IGr; */
404: row[0] = net_start + 2*gbus[i]+1;
405: col[0] = idx+2; col[1] = idx + 4; col[2] = idx + 5;
406: val[0] = -dIGr_ddelta; val[1] = -dIGr_dId; val[2] = -dIGr_dIq;
407: MatSetValues(J,1,row,3,col,val,INSERT_VALUES);
409: Vm = PetscSqrtScalar(Vd*Vd + Vq*Vq); Vm2 = Vm*Vm;
411: /* fgen[idx+6] = (KE[i]*Efd + SE - VR)/TE[i]; */
412: /* SE = k1[i]*PetscExpScalar(k2[i]*Efd); */
414: dSE_dEfd = k1[i]*k2[i]*PetscExpScalar(k2[i]*Efd);
416: row[0] = idx + 6;
417: col[0] = idx + 6; col[1] = idx + 8;
418: val[0] = (KE[i] + dSE_dEfd)/TE[i]; val[1] = -1/TE[i];
419: MatSetValues(J,1,row,2,col,val,INSERT_VALUES);
421: /* Exciter differential equations */
423: /* fgen[idx+7] = (RF - KF[i]*Efd/TF[i])/TF[i]; */
424: row[0] = idx + 7;
425: col[0] = idx + 6; col[1] = idx + 7;
426: val[0] = (-KF[i]/TF[i])/TF[i]; val[1] = 1/TF[i];
427: MatSetValues(J,1,row,2,col,val,INSERT_VALUES);
429: /* fgen[idx+8] = (VR - KA[i]*RF + KA[i]*KF[i]*Efd/TF[i] - KA[i]*(Vref[i] - Vm))/TA[i]; */
430: /* Vm = (Vd^2 + Vq^2)^0.5; */
432: dVm_dVd = Vd/Vm; dVm_dVq = Vq/Vm;
433: dVm_dVr = dVm_dVd*dVd_dVr + dVm_dVq*dVq_dVr;
434: dVm_dVi = dVm_dVd*dVd_dVi + dVm_dVq*dVq_dVi;
435: row[0] = idx + 8;
436: col[0] = idx + 6; col[1] = idx + 7; col[2] = idx + 8;
437: val[0] = (KA[i]*KF[i]/TF[i])/TA[i]; val[1] = -KA[i]/TA[i]; val[2] = 1/TA[i];
438: col[3] = net_start + 2*gbus[i]; col[4] = net_start + 2*gbus[i]+1;
439: val[3] = KA[i]*dVm_dVr/TA[i]; val[4] = KA[i]*dVm_dVi/TA[i];
440: MatSetValues(J,1,row,5,col,val,INSERT_VALUES);
441: idx = idx + 9;
442: }
444: for (i=0; i<nbus; i++) {
445: MatGetRow(user->Ybus,2*i,&ncols,&cols,&yvals);
446: row[0] = net_start + 2*i;
447: for (k=0; k<ncols; k++) {
448: col[k] = net_start + cols[k];
449: val[k] = yvals[k];
450: }
451: MatSetValues(J,1,row,ncols,col,val,INSERT_VALUES);
452: MatRestoreRow(user->Ybus,2*i,&ncols,&cols,&yvals);
454: MatGetRow(user->Ybus,2*i+1,&ncols,&cols,&yvals);
455: row[0] = net_start + 2*i+1;
456: for (k=0; k<ncols; k++) {
457: col[k] = net_start + cols[k];
458: val[k] = yvals[k];
459: }
460: MatSetValues(J,1,row,ncols,col,val,INSERT_VALUES);
461: MatRestoreRow(user->Ybus,2*i+1,&ncols,&cols,&yvals);
462: }
464: MatAssemblyBegin(J,MAT_FLUSH_ASSEMBLY);
465: MatAssemblyEnd(J,MAT_FLUSH_ASSEMBLY);
467: VecGetArray(user->V0,&v0);
468: for (i=0; i < nload; i++) {
469: Vr = xnet[2*lbus[i]]; /* Real part of load bus voltage */
470: Vi = xnet[2*lbus[i]+1]; /* Imaginary part of the load bus voltage */
471: Vm = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm; Vm4 = Vm2*Vm2;
472: Vm0 = PetscSqrtScalar(v0[2*lbus[i]]*v0[2*lbus[i]] + v0[2*lbus[i]+1]*v0[2*lbus[i]+1]);
473: PD = QD = 0.0;
474: dPD_dVr = dPD_dVi = dQD_dVr = dQD_dVi = 0.0;
475: for (k=0; k < ld_nsegsp[i]; k++) {
476: PD += ld_alphap[k]*PD0[i]*PetscPowScalar((Vm/Vm0),ld_betap[k]);
477: dPD_dVr += ld_alphap[k]*ld_betap[k]*PD0[i]*PetscPowScalar((1/Vm0),ld_betap[k])*Vr*PetscPowScalar(Vm,(ld_betap[k]-2));
478: dPD_dVi += ld_alphap[k]*ld_betap[k]*PD0[i]*PetscPowScalar((1/Vm0),ld_betap[k])*Vi*PetscPowScalar(Vm,(ld_betap[k]-2));
479: }
480: for (k=0; k < ld_nsegsq[i]; k++) {
481: QD += ld_alphaq[k]*QD0[i]*PetscPowScalar((Vm/Vm0),ld_betaq[k]);
482: dQD_dVr += ld_alphaq[k]*ld_betaq[k]*QD0[i]*PetscPowScalar((1/Vm0),ld_betaq[k])*Vr*PetscPowScalar(Vm,(ld_betaq[k]-2));
483: dQD_dVi += ld_alphaq[k]*ld_betaq[k]*QD0[i]*PetscPowScalar((1/Vm0),ld_betaq[k])*Vi*PetscPowScalar(Vm,(ld_betaq[k]-2));
484: }
486: /* IDr = (PD*Vr + QD*Vi)/Vm2; */
487: /* IDi = (-QD*Vr + PD*Vi)/Vm2; */
489: dIDr_dVr = (dPD_dVr*Vr + dQD_dVr*Vi + PD)/Vm2 - ((PD*Vr + QD*Vi)*2*Vr)/Vm4;
490: dIDr_dVi = (dPD_dVi*Vr + dQD_dVi*Vi + QD)/Vm2 - ((PD*Vr + QD*Vi)*2*Vi)/Vm4;
492: dIDi_dVr = (-dQD_dVr*Vr + dPD_dVr*Vi - QD)/Vm2 - ((-QD*Vr + PD*Vi)*2*Vr)/Vm4;
493: dIDi_dVi = (-dQD_dVi*Vr + dPD_dVi*Vi + PD)/Vm2 - ((-QD*Vr + PD*Vi)*2*Vi)/Vm4;
496: /* fnet[2*lbus[i]] += IDi; */
497: row[0] = net_start + 2*lbus[i];
498: col[0] = net_start + 2*lbus[i]; col[1] = net_start + 2*lbus[i]+1;
499: val[0] = dIDi_dVr; val[1] = dIDi_dVi;
500: MatSetValues(J,1,row,2,col,val,ADD_VALUES);
501: /* fnet[2*lbus[i]+1] += IDr; */
502: row[0] = net_start + 2*lbus[i]+1;
503: col[0] = net_start + 2*lbus[i]; col[1] = net_start + 2*lbus[i]+1;
504: val[0] = dIDr_dVr; val[1] = dIDr_dVi;
505: MatSetValues(J,1,row,2,col,val,ADD_VALUES);
506: }
507: VecRestoreArray(user->V0,&v0);
509: VecRestoreArray(Xgen,&xgen);
510: VecRestoreArray(Xnet,&xnet);
512: DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);
514: MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
515: MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
516: return(0);
517: }
519: int main(int argc,char **argv)
520: {
521: EPS eps;
522: EPSType type;
524: PetscMPIInt size;
525: Userctx user;
526: PetscViewer Xview,Ybusview;
527: Vec X,Xr,Xi;
528: Mat J,Jred=NULL;
529: IS is0,is1;
530: PetscInt i,*idx2,its,nev,nconv;
531: PetscReal error,re,im;
532: PetscScalar kr,ki;
533: PetscBool terse;
535: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
536: MPI_Comm_size(PETSC_COMM_WORLD,&size);
537: if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs");
538: /* show detailed info unless -terse option is given by user */
539: PetscOptionsHasName(NULL,NULL,"-terse",&terse);
541: user.neqs_gen = 9*ngen; /* # eqs. for generator subsystem */
542: user.neqs_net = 2*nbus; /* # eqs. for network subsystem */
543: user.neqs_pgrid = user.neqs_gen + user.neqs_net;
544: PetscPrintf(PETSC_COMM_WORLD,"\nStability analysis in a network with %D buses and %D generators\n\n",nbus,ngen);
546: /* Create indices for differential and algebraic equations */
547: PetscMalloc1(7*ngen,&idx2);
548: for (i=0; i<ngen; i++) {
549: 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;
550: idx2[7*i+4] = 9*i+6; idx2[7*i+5] = 9*i+7; idx2[7*i+6] = 9*i+8;
551: }
552: ISCreateGeneral(PETSC_COMM_WORLD,7*ngen,idx2,PETSC_COPY_VALUES,&user.is_diff);
553: ISComplement(user.is_diff,0,user.neqs_pgrid,&user.is_alg);
554: PetscFree(idx2);
556: /* Read initial voltage vector and Ybus */
557: PetscViewerBinaryOpen(PETSC_COMM_WORLD,"X.bin",FILE_MODE_READ,&Xview);
558: PetscViewerBinaryOpen(PETSC_COMM_WORLD,"Ybus.bin",FILE_MODE_READ,&Ybusview);
560: VecCreate(PETSC_COMM_WORLD,&user.V0);
561: VecSetSizes(user.V0,PETSC_DECIDE,user.neqs_net);
562: VecLoad(user.V0,Xview);
564: MatCreate(PETSC_COMM_WORLD,&user.Ybus);
565: MatSetSizes(user.Ybus,PETSC_DECIDE,PETSC_DECIDE,user.neqs_net,user.neqs_net);
566: MatSetType(user.Ybus,MATBAIJ);
567: /* MatSetBlockSize(user.Ybus,2); */
568: MatLoad(user.Ybus,Ybusview);
570: PetscViewerDestroy(&Xview);
571: PetscViewerDestroy(&Ybusview);
573: /* Create DMs for generator and network subsystems */
574: DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,user.neqs_gen,1,1,NULL,&user.dmgen);
575: DMSetOptionsPrefix(user.dmgen,"dmgen_");
576: DMSetFromOptions(user.dmgen);
577: DMSetUp(user.dmgen);
578: DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,user.neqs_net,1,1,NULL,&user.dmnet);
579: DMSetOptionsPrefix(user.dmnet,"dmnet_");
580: DMSetFromOptions(user.dmnet);
581: DMSetUp(user.dmnet);
583: /* Create a composite DM packer and add the two DMs */
584: DMCompositeCreate(PETSC_COMM_WORLD,&user.dmpgrid);
585: DMSetOptionsPrefix(user.dmpgrid,"pgrid_");
586: DMCompositeAddDM(user.dmpgrid,user.dmgen);
587: DMCompositeAddDM(user.dmpgrid,user.dmnet);
589: DMCreateGlobalVector(user.dmpgrid,&X);
591: MatCreate(PETSC_COMM_WORLD,&J);
592: MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,user.neqs_pgrid,user.neqs_pgrid);
593: MatSetFromOptions(J);
594: PreallocateJacobian(J,&user);
596: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
597: Set initial conditions
598: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
599: SetInitialGuess(X,&user);
601: /* Form Jacobian */
602: ResidualJacobian(X,J,(void*)&user);
603: MatScale(J,-1);
604: is0 = user.is_diff;
605: is1 = user.is_alg;
607: MatGetSchurComplement(J,is1,is1,is0,is0,MAT_IGNORE_MATRIX,NULL,MAT_SCHUR_COMPLEMENT_AINV_DIAG,MAT_INITIAL_MATRIX,&Jred);
609: if (!terse) {
610: MatView(Jred,NULL);
611: }
613: MatCreateVecs(Jred,NULL,&Xr);
614: MatCreateVecs(Jred,NULL,&Xi);
616: /* Create the eigensolver and set the various options */
617: EPSCreate(PETSC_COMM_WORLD,&eps);
618: EPSSetOperators(eps,Jred,NULL);
619: EPSSetProblemType(eps,EPS_NHEP);
620: EPSSetFromOptions(eps);
622: /* Solve the eigenvalue problem */
623: EPSSolve(eps);
625: EPSGetIterationNumber(eps,&its);
626: PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the eigensolver: %D\n",its);
627: EPSGetType(eps,&type);
628: PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n", type);
629: EPSGetDimensions(eps,&nev,NULL,NULL);
630: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);
632: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
633: Display solution and clean up
634: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
635: if (terse) {
636: EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
637: } else {
638: /* Get number of converged approximate eigenpairs */
639: EPSGetConverged(eps,&nconv);
640: PetscPrintf(PETSC_COMM_WORLD," Number of converged eigenpairs: %D\n\n",nconv);
642: if (nconv>0) {
643: /* Display eigenvalues and relative errors */
644: PetscPrintf(PETSC_COMM_WORLD,
645: " k ||Ax-kx||/||kx||\n"
646: " ----------------- ------------------\n");
648: for (i=0;i<nconv;i++) {
649: /* Get converged eigenpairs: i-th eigenvalue is stored in kr (real part) and
650: ki (imaginary part) */
651: EPSGetEigenpair(eps,i,&kr,&ki,Xr,Xi);
652: /* Compute the relative error associated to each eigenpair */
653: EPSComputeError(eps,i,EPS_ERROR_RELATIVE,&error);
655: #if defined(PETSC_USE_COMPLEX)
656: re = PetscRealPart(kr);
657: im = PetscImaginaryPart(kr);
658: #else
659: re = kr;
660: im = ki;
661: #endif
662: if (im!=0.0) {
663: PetscPrintf(PETSC_COMM_WORLD," %9f%+9fi %12g\n",(double)re,(double)im,(double)error);
664: } else {
665: PetscPrintf(PETSC_COMM_WORLD," %12f %12g\n",(double)re,(double)error);
666: }
667: }
668: PetscPrintf(PETSC_COMM_WORLD,"\n");
669: }
670: }
672: /* Free work space */
673: EPSDestroy(&eps);
674: MatDestroy(&J);
675: MatDestroy(&Jred);
676: MatDestroy(&user.Ybus);
677: VecDestroy(&X);
678: VecDestroy(&Xr);
679: VecDestroy(&Xi);
680: VecDestroy(&user.V0);
681: DMDestroy(&user.dmgen);
682: DMDestroy(&user.dmnet);
683: DMDestroy(&user.dmpgrid);
684: ISDestroy(&user.is_diff);
685: ISDestroy(&user.is_alg);
686: SlepcFinalize();
687: return ierr;
688: }
690: /*TEST
692: build:
693: requires: !complex
695: test:
696: suffix: 1
697: args: -terse
698: requires: double !complex !define(PETSC_USE_64BIT_INDICES)
699: localrunfiles: X.bin Ybus.bin
701: TEST*/