Actual source code: ex31.c

slepc-3.17.0 2022-03-31
Report Typos and Errors
  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:   *Fr =  Fd*PetscSinScalar(delta) + Fq*PetscCosScalar(delta);
140:   *Fi = -Fd*PetscCosScalar(delta) + Fq*PetscSinScalar(delta);
141:   PetscFunctionReturn(0);
142: }

144: /* Converts from network frame ([phase a real,imag) to machine (dq) reference frame */
145: PetscErrorCode ri2dq(PetscScalar Fr,PetscScalar Fi,PetscScalar delta,PetscScalar *Fd,PetscScalar *Fq)
146: {
147:   *Fd =  Fr*PetscSinScalar(delta) - Fi*PetscCosScalar(delta);
148:   *Fq =  Fr*PetscCosScalar(delta) + Fi*PetscSinScalar(delta);
149:   PetscFunctionReturn(0);
150: }

152: PetscErrorCode SetInitialGuess(Vec X,Userctx *user)
153: {
154:   Vec            Xgen,Xnet;
155:   PetscScalar    *xgen,*xnet;
156:   PetscInt       i,idx=0;
157:   PetscScalar    Vr,Vi,IGr,IGi,Vm,Vm2;
158:   PetscScalar    Eqp,Edp,delta;
159:   PetscScalar    Efd,RF,VR; /* Exciter variables */
160:   PetscScalar    Id,Iq;  /* Generator dq axis currents */
161:   PetscScalar    theta,Vd,Vq,SE;

163:   M[0] = 2*H[0]/w_s; M[1] = 2*H[1]/w_s; M[2] = 2*H[2]/w_s;
164:       /*      D[0] = 0.1*M[0]; D[1] = 0.1*M[1]; D[2] = 0.1*M[2];
165:        */
166:   D[0] = D[1] = D[2] = 0.0;
167:   DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);

169:   /* Network subsystem initialization */
170:   VecCopy(user->V0,Xnet);

172:   /* Generator subsystem initialization */
173:   VecGetArray(Xgen,&xgen);
174:   VecGetArray(Xnet,&xnet);

176:   for (i=0; i < ngen; i++) {
177:     Vr  = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */
178:     Vi  = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */
179:     Vm  = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm;
180:     IGr = (Vr*PG[i] + Vi*QG[i])/Vm2;
181:     IGi = (Vi*PG[i] - Vr*QG[i])/Vm2;

183:     delta = atan2(Vi+Xq[i]*IGr,Vr-Xq[i]*IGi); /* Machine angle */

185:     theta = PETSC_PI/2.0 - delta;

187:     Id = IGr*PetscCosScalar(theta) - IGi*PetscSinScalar(theta); /* d-axis stator current */
188:     Iq = IGr*PetscSinScalar(theta) + IGi*PetscCosScalar(theta); /* q-axis stator current */

190:     Vd = Vr*PetscCosScalar(theta) - Vi*PetscSinScalar(theta);
191:     Vq = Vr*PetscSinScalar(theta) + Vi*PetscCosScalar(theta);

193:     Edp = Vd + Rs[i]*Id - Xqp[i]*Iq; /* d-axis transient EMF */
194:     Eqp = Vq + Rs[i]*Iq + Xdp[i]*Id; /* q-axis transient EMF */

196:     TM[i] = PG[i];

198:     /* The generator variables are ordered as [Eqp,Edp,delta,w,Id,Iq] */
199:     xgen[idx]   = Eqp;
200:     xgen[idx+1] = Edp;
201:     xgen[idx+2] = delta;
202:     xgen[idx+3] = w_s;

204:     idx = idx + 4;

206:     xgen[idx]   = Id;
207:     xgen[idx+1] = Iq;

209:     idx = idx + 2;

211:     /* Exciter */
212:     Efd = Eqp + (Xd[i] - Xdp[i])*Id;
213:     SE  = k1[i]*PetscExpScalar(k2[i]*Efd);
214:     VR  =  KE[i]*Efd + SE;
215:     RF  =  KF[i]*Efd/TF[i];

217:     xgen[idx]   = Efd;
218:     xgen[idx+1] = RF;
219:     xgen[idx+2] = VR;

221:     Vref[i] = Vm + (VR/KA[i]);

223:     idx = idx + 3;
224:   }

226:   VecRestoreArray(Xgen,&xgen);
227:   VecRestoreArray(Xnet,&xnet);

229:   /* VecView(Xgen,0); */
230:   DMCompositeGather(user->dmpgrid,INSERT_VALUES,X,Xgen,Xnet);
231:   DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);
232:   PetscFunctionReturn(0);
233: }

235: PetscErrorCode PreallocateJacobian(Mat J,Userctx *user)
236: {
237:   PetscInt       *d_nnz;
238:   PetscInt       i,idx=0,start=0;
239:   PetscInt       ncols;

241:   PetscMalloc1(user->neqs_pgrid,&d_nnz);
242:   for (i=0; i<user->neqs_pgrid; i++) d_nnz[i] = 0;
243:   /* Generator subsystem */
244:   for (i=0; i < ngen; i++) {

246:     d_nnz[idx]   += 3;
247:     d_nnz[idx+1] += 2;
248:     d_nnz[idx+2] += 2;
249:     d_nnz[idx+3] += 5;
250:     d_nnz[idx+4] += 6;
251:     d_nnz[idx+5] += 6;

253:     d_nnz[user->neqs_gen+2*gbus[i]]   += 3;
254:     d_nnz[user->neqs_gen+2*gbus[i]+1] += 3;

256:     d_nnz[idx+6] += 2;
257:     d_nnz[idx+7] += 2;
258:     d_nnz[idx+8] += 5;

260:     idx = idx + 9;
261:   }

263:   start = user->neqs_gen;

265:   for (i=0; i < nbus; i++) {
266:     MatGetRow(user->Ybus,2*i,&ncols,NULL,NULL);
267:     d_nnz[start+2*i]   += ncols;
268:     d_nnz[start+2*i+1] += ncols;
269:     MatRestoreRow(user->Ybus,2*i,&ncols,NULL,NULL);
270:   }

272:   MatSeqAIJSetPreallocation(J,0,d_nnz);

274:   PetscFree(d_nnz);
275:   PetscFunctionReturn(0);
276: }

278: /*
279:    J = [-df_dx, -df_dy
280:         dg_dx, dg_dy]
281: */
282: PetscErrorCode ResidualJacobian(Vec X,Mat J,void *ctx)
283: {
284:   Userctx        *user=(Userctx*)ctx;
285:   Vec            Xgen,Xnet;
286:   PetscScalar    *xgen,*xnet;
287:   PetscInt       i,idx=0;
288:   PetscScalar    Vr,Vi,Vm,Vm2;
289:   PetscScalar    Eqp,Edp,delta; /* Generator variables */
290:   PetscScalar    Efd;
291:   PetscScalar    Id,Iq;  /* Generator dq axis currents */
292:   PetscScalar    Vd,Vq;
293:   PetscScalar    val[10];
294:   PetscInt       row[2],col[10];
295:   PetscInt       net_start=user->neqs_gen;
296:   PetscScalar    Zdq_inv[4],det;
297:   PetscScalar    dVd_dVr,dVd_dVi,dVq_dVr,dVq_dVi,dVd_ddelta,dVq_ddelta;
298:   PetscScalar    dIGr_ddelta,dIGi_ddelta,dIGr_dId,dIGr_dIq,dIGi_dId,dIGi_dIq;
299:   PetscScalar    dSE_dEfd;
300:   PetscScalar    dVm_dVd,dVm_dVq,dVm_dVr,dVm_dVi;
301:   PetscInt          ncols;
302:   const PetscInt    *cols;
303:   const PetscScalar *yvals;
304:   PetscInt          k;
305:   PetscScalar PD,QD,Vm0,*v0,Vm4;
306:   PetscScalar dPD_dVr,dPD_dVi,dQD_dVr,dQD_dVi;
307:   PetscScalar dIDr_dVr,dIDr_dVi,dIDi_dVr,dIDi_dVi;

309:   MatZeroEntries(J);
310:   DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
311:   DMCompositeScatter(user->dmpgrid,X,Xgen,Xnet);

313:   VecGetArray(Xgen,&xgen);
314:   VecGetArray(Xnet,&xnet);

316:   /* Generator subsystem */
317:   for (i=0; i < ngen; i++) {
318:     Eqp   = xgen[idx];
319:     Edp   = xgen[idx+1];
320:     delta = xgen[idx+2];
321:     Id    = xgen[idx+4];
322:     Iq    = xgen[idx+5];
323:     Efd   = xgen[idx+6];

325:     /*    fgen[idx]   = (Eqp + (Xd[i] - Xdp[i])*Id - Efd)/Td0p[i]; */
326:     row[0] = idx;
327:     col[0] = idx;           col[1] = idx+4;          col[2] = idx+6;
328:     val[0] = 1/ Td0p[i]; val[1] = (Xd[i] - Xdp[i])/ Td0p[i]; val[2] = -1/Td0p[i];

330:     MatSetValues(J,1,row,3,col,val,INSERT_VALUES);

332:     /*    fgen[idx+1] = (Edp - (Xq[i] - Xqp[i])*Iq)/Tq0p[i]; */
333:     row[0] = idx + 1;
334:     col[0] = idx + 1;       col[1] = idx+5;
335:     val[0] = 1/Tq0p[i]; val[1] = -(Xq[i] - Xqp[i])/Tq0p[i];
336:     MatSetValues(J,1,row,2,col,val,INSERT_VALUES);

338:     /*    fgen[idx+2] = - w + w_s; */
339:     row[0] = idx + 2;
340:     col[0] = idx + 2; col[1] = idx + 3;
341:     val[0] = 0;       val[1] = -1;
342:     MatSetValues(J,1,row,2,col,val,INSERT_VALUES);

344:     /*    fgen[idx+3] = (-TM[i] + Edp*Id + Eqp*Iq + (Xqp[i] - Xdp[i])*Id*Iq + D[i]*(w - w_s))/M[i]; */
345:     row[0] = idx + 3;
346:     col[0] = idx; col[1] = idx + 1; col[2] = idx + 3;       col[3] = idx + 4;                  col[4] = idx + 5;
347:     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];
348:     MatSetValues(J,1,row,5,col,val,INSERT_VALUES);

350:     Vr   = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */
351:     Vi   = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */
352:     ri2dq(Vr,Vi,delta,&Vd,&Vq);

354:     det = Rs[i]*Rs[i] + Xdp[i]*Xqp[i];

356:     Zdq_inv[0] = Rs[i]/det;
357:     Zdq_inv[1] = Xqp[i]/det;
358:     Zdq_inv[2] = -Xdp[i]/det;
359:     Zdq_inv[3] = Rs[i]/det;

361:     dVd_dVr    = PetscSinScalar(delta); dVd_dVi = -PetscCosScalar(delta);
362:     dVq_dVr    = PetscCosScalar(delta); dVq_dVi = PetscSinScalar(delta);
363:     dVd_ddelta = Vr*PetscCosScalar(delta) + Vi*PetscSinScalar(delta);
364:     dVq_ddelta = -Vr*PetscSinScalar(delta) + Vi*PetscCosScalar(delta);

366:     /*    fgen[idx+4] = Zdq_inv[0]*(-Edp + Vd) + Zdq_inv[1]*(-Eqp + Vq) + Id; */
367:     row[0] = idx+4;
368:     col[0] = idx;         col[1] = idx+1;        col[2] = idx + 2;
369:     val[0] = -Zdq_inv[1]; val[1] = -Zdq_inv[0];  val[2] = Zdq_inv[0]*dVd_ddelta + Zdq_inv[1]*dVq_ddelta;
370:     col[3] = idx + 4; col[4] = net_start+2*gbus[i];                     col[5] = net_start + 2*gbus[i]+1;
371:     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;
372:     MatSetValues(J,1,row,6,col,val,INSERT_VALUES);

374:     /*  fgen[idx+5] = Zdq_inv[2]*(-Edp + Vd) + Zdq_inv[3]*(-Eqp + Vq) + Iq; */
375:     row[0] = idx+5;
376:     col[0] = idx;         col[1] = idx+1;        col[2] = idx + 2;
377:     val[0] = -Zdq_inv[3]; val[1] = -Zdq_inv[2];  val[2] = Zdq_inv[2]*dVd_ddelta + Zdq_inv[3]*dVq_ddelta;
378:     col[3] = idx + 5; col[4] = net_start+2*gbus[i];                     col[5] = net_start + 2*gbus[i]+1;
379:     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;
380:     MatSetValues(J,1,row,6,col,val,INSERT_VALUES);

382:     dIGr_ddelta = Id*PetscCosScalar(delta) - Iq*PetscSinScalar(delta);
383:     dIGi_ddelta = Id*PetscSinScalar(delta) + Iq*PetscCosScalar(delta);
384:     dIGr_dId    = PetscSinScalar(delta);  dIGr_dIq = PetscCosScalar(delta);
385:     dIGi_dId    = -PetscCosScalar(delta); dIGi_dIq = PetscSinScalar(delta);

387:     /* fnet[2*gbus[i]]   -= IGi; */
388:     row[0] = net_start + 2*gbus[i];
389:     col[0] = idx+2;        col[1] = idx + 4;   col[2] = idx + 5;
390:     val[0] = -dIGi_ddelta; val[1] = -dIGi_dId; val[2] = -dIGi_dIq;
391:     MatSetValues(J,1,row,3,col,val,INSERT_VALUES);

393:     /* fnet[2*gbus[i]+1]   -= IGr; */
394:     row[0] = net_start + 2*gbus[i]+1;
395:     col[0] = idx+2;        col[1] = idx + 4;   col[2] = idx + 5;
396:     val[0] = -dIGr_ddelta; val[1] = -dIGr_dId; val[2] = -dIGr_dIq;
397:     MatSetValues(J,1,row,3,col,val,INSERT_VALUES);

399:     Vm = PetscSqrtScalar(Vd*Vd + Vq*Vq); Vm2 = Vm*Vm;

401:     /*    fgen[idx+6] = (KE[i]*Efd + SE - VR)/TE[i]; */
402:     /*    SE  = k1[i]*PetscExpScalar(k2[i]*Efd); */

404:     dSE_dEfd = k1[i]*k2[i]*PetscExpScalar(k2[i]*Efd);

406:     row[0] = idx + 6;
407:     col[0] = idx + 6;                     col[1] = idx + 8;
408:     val[0] = (KE[i] + dSE_dEfd)/TE[i];  val[1] = -1/TE[i];
409:     MatSetValues(J,1,row,2,col,val,INSERT_VALUES);

411:     /* Exciter differential equations */

413:     /*    fgen[idx+7] = (RF - KF[i]*Efd/TF[i])/TF[i]; */
414:     row[0] = idx + 7;
415:     col[0] = idx + 6;       col[1] = idx + 7;
416:     val[0] = (-KF[i]/TF[i])/TF[i];  val[1] = 1/TF[i];
417:     MatSetValues(J,1,row,2,col,val,INSERT_VALUES);

419:     /*    fgen[idx+8] = (VR - KA[i]*RF + KA[i]*KF[i]*Efd/TF[i] - KA[i]*(Vref[i] - Vm))/TA[i]; */
420:     /* Vm = (Vd^2 + Vq^2)^0.5; */

422:     dVm_dVd = Vd/Vm; dVm_dVq = Vq/Vm;
423:     dVm_dVr = dVm_dVd*dVd_dVr + dVm_dVq*dVq_dVr;
424:     dVm_dVi = dVm_dVd*dVd_dVi + dVm_dVq*dVq_dVi;
425:     row[0]  = idx + 8;
426:     col[0]  = idx + 6;           col[1] = idx + 7; col[2] = idx + 8;
427:     val[0]  = (KA[i]*KF[i]/TF[i])/TA[i]; val[1] = -KA[i]/TA[i];  val[2] = 1/TA[i];
428:     col[3]  = net_start + 2*gbus[i]; col[4] = net_start + 2*gbus[i]+1;
429:     val[3]  = KA[i]*dVm_dVr/TA[i];         val[4] = KA[i]*dVm_dVi/TA[i];
430:     MatSetValues(J,1,row,5,col,val,INSERT_VALUES);
431:     idx     = idx + 9;
432:   }

434:   for (i=0; i<nbus; i++) {
435:     MatGetRow(user->Ybus,2*i,&ncols,&cols,&yvals);
436:     row[0] = net_start + 2*i;
437:     for (k=0; k<ncols; k++) {
438:       col[k] = net_start + cols[k];
439:       val[k] = yvals[k];
440:     }
441:     MatSetValues(J,1,row,ncols,col,val,INSERT_VALUES);
442:     MatRestoreRow(user->Ybus,2*i,&ncols,&cols,&yvals);

444:     MatGetRow(user->Ybus,2*i+1,&ncols,&cols,&yvals);
445:     row[0] = net_start + 2*i+1;
446:     for (k=0; k<ncols; k++) {
447:       col[k] = net_start + cols[k];
448:       val[k] = yvals[k];
449:     }
450:     MatSetValues(J,1,row,ncols,col,val,INSERT_VALUES);
451:     MatRestoreRow(user->Ybus,2*i+1,&ncols,&cols,&yvals);
452:   }

454:   MatAssemblyBegin(J,MAT_FLUSH_ASSEMBLY);
455:   MatAssemblyEnd(J,MAT_FLUSH_ASSEMBLY);

457:   VecGetArray(user->V0,&v0);
458:   for (i=0; i < nload; i++) {
459:     Vr      = xnet[2*lbus[i]]; /* Real part of load bus voltage */
460:     Vi      = xnet[2*lbus[i]+1]; /* Imaginary part of the load bus voltage */
461:     Vm      = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm; Vm4 = Vm2*Vm2;
462:     Vm0     = PetscSqrtScalar(v0[2*lbus[i]]*v0[2*lbus[i]] + v0[2*lbus[i]+1]*v0[2*lbus[i]+1]);
463:     PD      = QD = 0.0;
464:     dPD_dVr = dPD_dVi = dQD_dVr = dQD_dVi = 0.0;
465:     for (k=0; k < ld_nsegsp[i]; k++) {
466:       PD      += ld_alphap[k]*PD0[i]*PetscPowScalar((Vm/Vm0),ld_betap[k]);
467:       dPD_dVr += ld_alphap[k]*ld_betap[k]*PD0[i]*PetscPowScalar((1/Vm0),ld_betap[k])*Vr*PetscPowScalar(Vm,(ld_betap[k]-2));
468:       dPD_dVi += ld_alphap[k]*ld_betap[k]*PD0[i]*PetscPowScalar((1/Vm0),ld_betap[k])*Vi*PetscPowScalar(Vm,(ld_betap[k]-2));
469:     }
470:     for (k=0; k < ld_nsegsq[i]; k++) {
471:       QD      += ld_alphaq[k]*QD0[i]*PetscPowScalar((Vm/Vm0),ld_betaq[k]);
472:       dQD_dVr += ld_alphaq[k]*ld_betaq[k]*QD0[i]*PetscPowScalar((1/Vm0),ld_betaq[k])*Vr*PetscPowScalar(Vm,(ld_betaq[k]-2));
473:       dQD_dVi += ld_alphaq[k]*ld_betaq[k]*QD0[i]*PetscPowScalar((1/Vm0),ld_betaq[k])*Vi*PetscPowScalar(Vm,(ld_betaq[k]-2));
474:     }

476:     /*    IDr = (PD*Vr + QD*Vi)/Vm2; */
477:     /*    IDi = (-QD*Vr + PD*Vi)/Vm2; */

479:     dIDr_dVr = (dPD_dVr*Vr + dQD_dVr*Vi + PD)/Vm2 - ((PD*Vr + QD*Vi)*2*Vr)/Vm4;
480:     dIDr_dVi = (dPD_dVi*Vr + dQD_dVi*Vi + QD)/Vm2 - ((PD*Vr + QD*Vi)*2*Vi)/Vm4;

482:     dIDi_dVr = (-dQD_dVr*Vr + dPD_dVr*Vi - QD)/Vm2 - ((-QD*Vr + PD*Vi)*2*Vr)/Vm4;
483:     dIDi_dVi = (-dQD_dVi*Vr + dPD_dVi*Vi + PD)/Vm2 - ((-QD*Vr + PD*Vi)*2*Vi)/Vm4;

485:     /*    fnet[2*lbus[i]]   += IDi; */
486:     row[0] = net_start + 2*lbus[i];
487:     col[0] = net_start + 2*lbus[i];  col[1] = net_start + 2*lbus[i]+1;
488:     val[0] = dIDi_dVr;               val[1] = dIDi_dVi;
489:     MatSetValues(J,1,row,2,col,val,ADD_VALUES);
490:     /*    fnet[2*lbus[i]+1] += IDr; */
491:     row[0] = net_start + 2*lbus[i]+1;
492:     col[0] = net_start + 2*lbus[i];  col[1] = net_start + 2*lbus[i]+1;
493:     val[0] = dIDr_dVr;               val[1] = dIDr_dVi;
494:     MatSetValues(J,1,row,2,col,val,ADD_VALUES);
495:   }
496:   VecRestoreArray(user->V0,&v0);

498:   VecRestoreArray(Xgen,&xgen);
499:   VecRestoreArray(Xnet,&xnet);

501:   DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);

503:   MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
504:   MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
505:   PetscFunctionReturn(0);
506: }

508: int main(int argc,char **argv)
509: {
510:   EPS            eps;
511:   EPSType        type;
512:   PetscMPIInt    size;
513:   Userctx        user;
514:   PetscViewer    Xview,Ybusview;
515:   Vec            X,Xr,Xi;
516:   Mat            J,Jred=NULL;
517:   IS             is0,is1;
518:   PetscInt       i,*idx2,its,nev,nconv;
519:   PetscReal      error,re,im;
520:   PetscScalar    kr,ki;
521:   PetscBool      terse;

523:   SlepcInitialize(&argc,&argv,(char*)0,help);
524:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
526:   /* show detailed info unless -terse option is given by user */
527:   PetscOptionsHasName(NULL,NULL,"-terse",&terse);

529:   user.neqs_gen   = 9*ngen; /* # eqs. for generator subsystem */
530:   user.neqs_net   = 2*nbus; /* # eqs. for network subsystem   */
531:   user.neqs_pgrid = user.neqs_gen + user.neqs_net;
532:   PetscPrintf(PETSC_COMM_WORLD,"\nStability analysis in a network with %" PetscInt_FMT " buses and %" PetscInt_FMT " generators\n\n",nbus,ngen);

534:   /* Create indices for differential and algebraic equations */
535:   PetscMalloc1(7*ngen,&idx2);
536:   for (i=0; i<ngen; i++) {
537:     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;
538:     idx2[7*i+4] = 9*i+6; idx2[7*i+5] = 9*i+7; idx2[7*i+6] = 9*i+8;
539:   }
540:   ISCreateGeneral(PETSC_COMM_WORLD,7*ngen,idx2,PETSC_COPY_VALUES,&user.is_diff);
541:   ISComplement(user.is_diff,0,user.neqs_pgrid,&user.is_alg);
542:   PetscFree(idx2);

544:   /* Read initial voltage vector and Ybus */
545:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"X.bin",FILE_MODE_READ,&Xview);
546:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"Ybus.bin",FILE_MODE_READ,&Ybusview);

548:   VecCreate(PETSC_COMM_WORLD,&user.V0);
549:   VecSetSizes(user.V0,PETSC_DECIDE,user.neqs_net);
550:   VecLoad(user.V0,Xview);

552:   MatCreate(PETSC_COMM_WORLD,&user.Ybus);
553:   MatSetSizes(user.Ybus,PETSC_DECIDE,PETSC_DECIDE,user.neqs_net,user.neqs_net);
554:   MatSetType(user.Ybus,MATBAIJ);
555:   /*  MatSetBlockSize(user.Ybus,2); */
556:   MatLoad(user.Ybus,Ybusview);

558:   PetscViewerDestroy(&Xview);
559:   PetscViewerDestroy(&Ybusview);

561:   /* Create DMs for generator and network subsystems */
562:   DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,user.neqs_gen,1,1,NULL,&user.dmgen);
563:   DMSetOptionsPrefix(user.dmgen,"dmgen_");
564:   DMSetFromOptions(user.dmgen);
565:   DMSetUp(user.dmgen);
566:   DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,user.neqs_net,1,1,NULL,&user.dmnet);
567:   DMSetOptionsPrefix(user.dmnet,"dmnet_");
568:   DMSetFromOptions(user.dmnet);
569:   DMSetUp(user.dmnet);

571:   /* Create a composite DM packer and add the two DMs */
572:   DMCompositeCreate(PETSC_COMM_WORLD,&user.dmpgrid);
573:   DMSetOptionsPrefix(user.dmpgrid,"pgrid_");
574:   DMCompositeAddDM(user.dmpgrid,user.dmgen);
575:   DMCompositeAddDM(user.dmpgrid,user.dmnet);

577:   DMCreateGlobalVector(user.dmpgrid,&X);

579:   MatCreate(PETSC_COMM_WORLD,&J);
580:   MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,user.neqs_pgrid,user.neqs_pgrid);
581:   MatSetFromOptions(J);
582:   PreallocateJacobian(J,&user);

584:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
585:      Set initial conditions
586:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
587:   SetInitialGuess(X,&user);

589:   /* Form Jacobian */
590:   ResidualJacobian(X,J,(void*)&user);
591:   MatScale(J,-1);
592:   is0 = user.is_diff;
593:   is1 = user.is_alg;

595:   MatGetSchurComplement(J,is1,is1,is0,is0,MAT_IGNORE_MATRIX,NULL,MAT_SCHUR_COMPLEMENT_AINV_DIAG,MAT_INITIAL_MATRIX,&Jred);

597:   if (!terse) MatView(Jred,NULL);

599:   MatCreateVecs(Jred,NULL,&Xr);
600:   MatCreateVecs(Jred,NULL,&Xi);

602:   /* Create the eigensolver and set the various options */
603:   EPSCreate(PETSC_COMM_WORLD,&eps);
604:   EPSSetOperators(eps,Jred,NULL);
605:   EPSSetProblemType(eps,EPS_NHEP);
606:   EPSSetFromOptions(eps);

608:   /* Solve the eigenvalue problem */
609:   EPSSolve(eps);

611:   EPSGetIterationNumber(eps,&its);
612:   PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the eigensolver: %" PetscInt_FMT "\n",its);
613:   EPSGetType(eps,&type);
614:   PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n", type);
615:   EPSGetDimensions(eps,&nev,NULL,NULL);
616:   PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev);

618:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
619:                     Display solution and clean up
620:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
621:   if (terse) EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
622:   else {
623:     /* Get number of converged approximate eigenpairs */
624:     EPSGetConverged(eps,&nconv);
625:     PetscPrintf(PETSC_COMM_WORLD," Number of converged eigenpairs: %" PetscInt_FMT "\n\n",nconv);

627:     if (nconv>0) {
628:       /* Display eigenvalues and relative errors */
629:       PetscCall(PetscPrintf(PETSC_COMM_WORLD,
630:            "           k          ||Ax-kx||/||kx||\n"
631:            "   ----------------- ------------------\n"));

633:       for (i=0;i<nconv;i++) {
634:         /* Get converged eigenpairs: i-th eigenvalue is stored in kr (real part) and
635:           ki (imaginary part) */
636:         EPSGetEigenpair(eps,i,&kr,&ki,Xr,Xi);
637:         /* Compute the relative error associated to each eigenpair */
638:         EPSComputeError(eps,i,EPS_ERROR_RELATIVE,&error);

640: #if defined(PETSC_USE_COMPLEX)
641:         re = PetscRealPart(kr);
642:         im = PetscImaginaryPart(kr);
643: #else
644:         re = kr;
645:         im = ki;
646: #endif
647:         if (im!=0.0) PetscPrintf(PETSC_COMM_WORLD," %9f%+9fi %12g\n",(double)re,(double)im,(double)error);
648:         else PetscPrintf(PETSC_COMM_WORLD,"   %12f       %12g\n",(double)re,(double)error);
649:       }
650:       PetscPrintf(PETSC_COMM_WORLD,"\n");
651:     }
652:   }

654:   /* Free work space */
655:   EPSDestroy(&eps);
656:   MatDestroy(&J);
657:   MatDestroy(&Jred);
658:   MatDestroy(&user.Ybus);
659:   VecDestroy(&X);
660:   VecDestroy(&Xr);
661:   VecDestroy(&Xi);
662:   VecDestroy(&user.V0);
663:   DMDestroy(&user.dmgen);
664:   DMDestroy(&user.dmnet);
665:   DMDestroy(&user.dmpgrid);
666:   ISDestroy(&user.is_diff);
667:   ISDestroy(&user.is_alg);
668:   SlepcFinalize();
669:   return 0;
670: }

672: /*TEST

674:    build:
675:       requires: !complex

677:    test:
678:       suffix: 1
679:       args: -terse
680:       requires: double !complex !defined(PETSC_USE_64BIT_INDICES)
681:       localrunfiles: X.bin Ybus.bin

683: TEST*/