Actual source code: ex31.c

slepc-3.7.0 2016-05-16
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.

  8:    SLEPc is free software: you can redistribute it and/or modify it under  the
  9:    terms of version 3 of the GNU Lesser General Public License as published by
 10:    the Free Software Foundation.

 12:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
 13:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
 14:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
 15:    more details.

 17:    You  should have received a copy of the GNU Lesser General  Public  License
 18:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 19:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 20: */

 22: static char help[] = "Power grid small signal stability analysis of WECC 9 bus system.\n\
 23: This example is based on the 9-bus (node) example given in the book Power\n\
 24: Systems Dynamics and Stability (Chapter 8) by P. Sauer and M. A. Pai.\n\
 25: The power grid in this example consists of 9 buses (nodes), 3 generators,\n\
 26: 3 loads, and 9 transmission lines. The network equations are written\n\
 27: in current balance form using rectangular coordinates. It uses the SLEPc\n\
 28: package to calculate the eigenvalues for small signal stability analysis\n\n";

 30: /*
 31:    This example is based on PETSc's ex9bus example (under TS).

 33:    The equations for the stability analysis are described by the DAE

 35:    \dot{x} = f(x,y,t)
 36:      0     = g(x,y,t)

 38:    where the generators are described by differential equations, while the algebraic
 39:    constraints define the network equations.

 41:    The generators are modeled with a 4th order differential equation describing the electrical
 42:    and mechanical dynamics. Each generator also has an exciter system modeled by 3rd order
 43:    diff. eqns. describing the exciter, voltage regulator, and the feedback stabilizer
 44:    mechanism.

 46:    The network equations are described by nodal current balance equations.
 47:     I(x,y) - Y*V = 0

 49:    where:
 50:     I(x,y) is the current injected from generators and loads.
 51:       Y    is the admittance matrix, and
 52:       V    is the voltage vector

 54:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 55:  
 56:    The linearized equations for the eigenvalue analysis are

 58:      \dot{\delta{x}} = f_x*\delta{x} + f_y*\delta{y}
 59:              0       = g_x*\delta{x} + g_y*\delta{y}

 61:    This gives the linearized sensitivity matrix
 62:      A = | f_x  f_y |
 63:          | g_x  g_y |

 65:    We are interested in the eigenvalues of the Schur complement of A
 66:      \hat{A} = f_x - g_x*inv(g_y)*f_y


 69:    Example contributed by: Shrirang Abhyankar
 70: */

 72: #include <petscdm.h>
 73: #include <petscdmda.h>
 74: #include <petscdmcomposite.h>
 75: #include <slepceps.h>

 77: #define freq 60
 78: #define w_s (2*PETSC_PI*freq)

 80: /* Sizes and indices */
 81: const PetscInt nbus    = 9; /* Number of network buses */
 82: const PetscInt ngen    = 3; /* Number of generators */
 83: const PetscInt nload   = 3; /* Number of loads */
 84: const PetscInt gbus[3] = {0,1,2}; /* Buses at which generators are incident */
 85: const PetscInt lbus[3] = {4,5,7}; /* Buses at which loads are incident */

 87: /* Generator real and reactive powers (found via loadflow) */
 88: const PetscScalar PG[3] = {0.716786142395021,1.630000000000000,0.850000000000000};
 89: const PetscScalar QG[3] = {0.270702180178785,0.066120127797275,-0.108402221791588};
 90: /* Generator constants */
 91: const PetscScalar H[3]    = {23.64,6.4,3.01};   /* Inertia constant */
 92: const PetscScalar Rs[3]   = {0.0,0.0,0.0}; /* Stator Resistance */
 93: const PetscScalar Xd[3]   = {0.146,0.8958,1.3125};  /* d-axis reactance */
 94: const PetscScalar Xdp[3]  = {0.0608,0.1198,0.1813}; /* d-axis transient reactance */
 95: 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 */
 96: const PetscScalar Xqp[3]  = {0.0969,0.1969,0.25}; /* q-axis transient reactance */
 97: const PetscScalar Td0p[3] = {8.96,6.0,5.89}; /* d-axis open circuit time constant */
 98: const PetscScalar Tq0p[3] = {0.31,0.535,0.6}; /* q-axis open circuit time constant */
 99: PetscScalar M[3]; /* M = 2*H/w_s */
100: PetscScalar D[3]; /* D = 0.1*M */

102: PetscScalar TM[3]; /* Mechanical Torque */
103: /* Exciter system constants */
104: const PetscScalar KA[3] = {20.0,20.0,20.0};  /* Voltage regulartor gain constant */
105: const PetscScalar TA[3] = {0.2,0.2,0.2};     /* Voltage regulator time constant */
106: const PetscScalar KE[3] = {1.0,1.0,1.0};     /* Exciter gain constant */
107: const PetscScalar TE[3] = {0.314,0.314,0.314}; /* Exciter time constant */
108: const PetscScalar KF[3] = {0.063,0.063,0.063};  /* Feedback stabilizer gain constant */
109: const PetscScalar TF[3] = {0.35,0.35,0.35};    /* Feedback stabilizer time constant */
110: const PetscScalar k1[3] = {0.0039,0.0039,0.0039};
111: const PetscScalar k2[3] = {1.555,1.555,1.555};  /* k1 and k2 for calculating the saturation function SE = k1*exp(k2*Efd) */

113: PetscScalar Vref[3];
114: /* Load constants
115:   We use a composite load model that describes the load and reactive powers at each time instant as follows
116:   P(t) = \sum\limits_{i=0}^ld_nsegsp \ld_alphap_i*P_D0(\frac{V_m(t)}{V_m0})^\ld_betap_i
117:   Q(t) = \sum\limits_{i=0}^ld_nsegsq \ld_alphaq_i*Q_D0(\frac{V_m(t)}{V_m0})^\ld_betaq_i
118:   where
119:     ld_nsegsp,ld_nsegsq - Number of individual load models for real and reactive power loads
120:     ld_alphap,ld_alphap - Percentage contribution (weights) or loads
121:     P_D0                - Real power load
122:     Q_D0                - Reactive power load
123:     V_m(t)              - Voltage magnitude at time t
124:     V_m0                - Voltage magnitude at t = 0
125:     ld_betap, ld_betaq  - exponents describing the load model for real and reactive part

127:     Note: All loads have the same characteristic currently.
128: */
129: const PetscScalar PD0[3] = {1.25,0.9,1.0};
130: const PetscScalar QD0[3] = {0.5,0.3,0.35};
131: const PetscInt    ld_nsegsp[3] = {3,3,3};
132: const PetscScalar ld_alphap[3] = {0.0,0.0,1.0};
133: const PetscScalar ld_betap[3]  = {2.0,1.0,0.0};
134: const PetscInt    ld_nsegsq[3] = {3,3,3};
135: const PetscScalar ld_alphaq[3] = {0.0,0.0,1.0};
136: const PetscScalar ld_betaq[3]  = {2.0,1.0,0.0};

138: typedef struct {
139:   DM       dmgen, dmnet; /* DMs to manage generator and network subsystem */
140:   DM       dmpgrid;      /* Composite DM to manage the entire power grid */
141:   Mat      Ybus;         /* Network admittance matrix */
142:   Vec      V0;           /* Initial voltage vector (Power flow solution) */
143:   PetscInt neqs_gen,neqs_net,neqs_pgrid;
144:   IS       is_diff;      /* indices for differential equations */
145:   IS       is_alg;       /* indices for algebraic equations */
146: } Userctx;

148: /* Converts from machine frame (dq) to network (phase a real,imag) reference frame */
151: PetscErrorCode dq2ri(PetscScalar Fd,PetscScalar Fq,PetscScalar delta,PetscScalar *Fr,PetscScalar *Fi)
152: {
154:   *Fr =  Fd*PetscSinScalar(delta) + Fq*PetscCosScalar(delta);
155:   *Fi = -Fd*PetscCosScalar(delta) + Fq*PetscSinScalar(delta);
156:   return(0);
157: }

159: /* Converts from network frame ([phase a real,imag) to machine (dq) reference frame */
162: PetscErrorCode ri2dq(PetscScalar Fr,PetscScalar Fi,PetscScalar delta,PetscScalar *Fd,PetscScalar *Fq)
163: {
165:   *Fd =  Fr*PetscSinScalar(delta) - Fi*PetscCosScalar(delta);
166:   *Fq =  Fr*PetscCosScalar(delta) + Fi*PetscSinScalar(delta);
167:   return(0);
168: }

172: PetscErrorCode SetInitialGuess(Vec X,Userctx *user)
173: {
175:   Vec            Xgen,Xnet;
176:   PetscScalar    *xgen,*xnet;
177:   PetscInt       i,idx=0;
178:   PetscScalar    Vr,Vi,IGr,IGi,Vm,Vm2;
179:   PetscScalar    Eqp,Edp,delta;
180:   PetscScalar    Efd,RF,VR; /* Exciter variables */
181:   PetscScalar    Id,Iq;  /* Generator dq axis currents */
182:   PetscScalar    theta,Vd,Vq,SE;

185:   M[0] = 2*H[0]/w_s; M[1] = 2*H[1]/w_s; M[2] = 2*H[2]/w_s;
186:       /*      D[0] = 0.1*M[0]; D[1] = 0.1*M[1]; D[2] = 0.1*M[2];
187:        */
188:   D[0] = D[1] = D[2] = 0.0;
189:   DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);

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

194:   /* Generator subsystem initialization */
195:   VecGetArray(Xgen,&xgen);
196:   VecGetArray(Xnet,&xnet);

198:   for (i=0; i < ngen; i++) {
199:     Vr  = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */
200:     Vi  = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */
201:     Vm  = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm;
202:     IGr = (Vr*PG[i] + Vi*QG[i])/Vm2;
203:     IGi = (Vi*PG[i] - Vr*QG[i])/Vm2;

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

207:     theta = PETSC_PI/2.0 - delta;

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

212:     Vd = Vr*PetscCosScalar(theta) - Vi*PetscSinScalar(theta);
213:     Vq = Vr*PetscSinScalar(theta) + Vi*PetscCosScalar(theta);

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

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

220:     /* The generator variables are ordered as [Eqp,Edp,delta,w,Id,Iq] */
221:     xgen[idx]   = Eqp;
222:     xgen[idx+1] = Edp;
223:     xgen[idx+2] = delta;
224:     xgen[idx+3] = w_s;

226:     idx = idx + 4;

228:     xgen[idx]   = Id;
229:     xgen[idx+1] = Iq;

231:     idx = idx + 2;

233:     /* Exciter */
234:     Efd = Eqp + (Xd[i] - Xdp[i])*Id;
235:     SE  = k1[i]*PetscExpScalar(k2[i]*Efd);
236:     VR  =  KE[i]*Efd + SE;
237:     RF  =  KF[i]*Efd/TF[i];

239:     xgen[idx]   = Efd;
240:     xgen[idx+1] = RF;
241:     xgen[idx+2] = VR;

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

245:     idx = idx + 3;
246:   }

248:   VecRestoreArray(Xgen,&xgen);
249:   VecRestoreArray(Xnet,&xnet);

251:   /* VecView(Xgen,0); */
252:   DMCompositeGather(user->dmpgrid,X,INSERT_VALUES,Xgen,Xnet);
253:   DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);
254:   return(0);
255: }

259: PetscErrorCode PreallocateJacobian(Mat J,Userctx *user)
260: {
262:   PetscInt       *d_nnz;
263:   PetscInt       i,idx=0,start=0;
264:   PetscInt       ncols;

267:   PetscMalloc1(user->neqs_pgrid,&d_nnz);
268:   for (i=0; i<user->neqs_pgrid; i++) d_nnz[i] = 0;
269:   /* Generator subsystem */
270:   for (i=0; i < ngen; i++) {

272:     d_nnz[idx]   += 3;
273:     d_nnz[idx+1] += 2;
274:     d_nnz[idx+2] += 2;
275:     d_nnz[idx+3] += 5;
276:     d_nnz[idx+4] += 6;
277:     d_nnz[idx+5] += 6;

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

282:     d_nnz[idx+6] += 2;
283:     d_nnz[idx+7] += 2;
284:     d_nnz[idx+8] += 5;

286:     idx = idx + 9;
287:   }

289:   start = user->neqs_gen;

291:   for (i=0; i < nbus; i++) {
292:     MatGetRow(user->Ybus,2*i,&ncols,NULL,NULL);
293:     d_nnz[start+2*i]   += ncols;
294:     d_nnz[start+2*i+1] += ncols;
295:     MatRestoreRow(user->Ybus,2*i,&ncols,NULL,NULL);
296:   }

298:   MatSeqAIJSetPreallocation(J,0,d_nnz);

300:   PetscFree(d_nnz);
301:   return(0);
302: }

304: /*
305:    J = [-df_dx, -df_dy
306:         dg_dx, dg_dy]
307: */
310: PetscErrorCode ResidualJacobian(Vec X,Mat J,void *ctx)
311: {
313:   Userctx        *user=(Userctx*)ctx;
314:   Vec            Xgen,Xnet;
315:   PetscScalar    *xgen,*xnet;
316:   PetscInt       i,idx=0;
317:   PetscScalar    Vr,Vi,Vm,Vm2;
318:   PetscScalar    Eqp,Edp,delta; /* Generator variables */
319:   PetscScalar    Efd;
320:   PetscScalar    Id,Iq;  /* Generator dq axis currents */
321:   PetscScalar    Vd,Vq;
322:   PetscScalar    val[10];
323:   PetscInt       row[2],col[10];
324:   PetscInt       net_start=user->neqs_gen;
325:   PetscScalar    Zdq_inv[4],det;
326:   PetscScalar    dVd_dVr,dVd_dVi,dVq_dVr,dVq_dVi,dVd_ddelta,dVq_ddelta;
327:   PetscScalar    dIGr_ddelta,dIGi_ddelta,dIGr_dId,dIGr_dIq,dIGi_dId,dIGi_dIq;
328:   PetscScalar    dSE_dEfd;
329:   PetscScalar    dVm_dVd,dVm_dVq,dVm_dVr,dVm_dVi;
330:   PetscInt          ncols;
331:   const PetscInt    *cols;
332:   const PetscScalar *yvals;
333:   PetscInt          k;
334:   PetscScalar PD,QD,Vm0,*v0,Vm4;
335:   PetscScalar dPD_dVr,dPD_dVi,dQD_dVr,dQD_dVi;
336:   PetscScalar dIDr_dVr,dIDr_dVi,dIDi_dVr,dIDi_dVi;


340:   MatZeroEntries(J);
341:   DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
342:   DMCompositeScatter(user->dmpgrid,X,Xgen,Xnet);

344:   VecGetArray(Xgen,&xgen);
345:   VecGetArray(Xnet,&xnet);

347:   /* Generator subsystem */
348:   for (i=0; i < ngen; i++) {
349:     Eqp   = xgen[idx];
350:     Edp   = xgen[idx+1];
351:     delta = xgen[idx+2];
352:     Id    = xgen[idx+4];
353:     Iq    = xgen[idx+5];
354:     Efd   = xgen[idx+6];

356:     /*    fgen[idx]   = (Eqp + (Xd[i] - Xdp[i])*Id - Efd)/Td0p[i]; */
357:     row[0] = idx;
358:     col[0] = idx;           col[1] = idx+4;          col[2] = idx+6;
359:     val[0] = 1/ Td0p[i]; val[1] = (Xd[i] - Xdp[i])/ Td0p[i]; val[2] = -1/Td0p[i];

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

363:     /*    fgen[idx+1] = (Edp - (Xq[i] - Xqp[i])*Iq)/Tq0p[i]; */
364:     row[0] = idx + 1;
365:     col[0] = idx + 1;       col[1] = idx+5;
366:     val[0] = 1/Tq0p[i]; val[1] = -(Xq[i] - Xqp[i])/Tq0p[i];
367:     MatSetValues(J,1,row,2,col,val,INSERT_VALUES);

369:     /*    fgen[idx+2] = - w + w_s; */
370:     row[0] = idx + 2;
371:     col[0] = idx + 2; col[1] = idx + 3;
372:     val[0] = 0;       val[1] = -1;
373:     MatSetValues(J,1,row,2,col,val,INSERT_VALUES);

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

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

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

387:     Zdq_inv[0] = Rs[i]/det;
388:     Zdq_inv[1] = Xqp[i]/det;
389:     Zdq_inv[2] = -Xdp[i]/det;
390:     Zdq_inv[3] = Rs[i]/det;

392:     dVd_dVr    = PetscSinScalar(delta); dVd_dVi = -PetscCosScalar(delta);
393:     dVq_dVr    = PetscCosScalar(delta); dVq_dVi = PetscSinScalar(delta);
394:     dVd_ddelta = Vr*PetscCosScalar(delta) + Vi*PetscSinScalar(delta);
395:     dVq_ddelta = -Vr*PetscSinScalar(delta) + Vi*PetscCosScalar(delta);

397:     /*    fgen[idx+4] = Zdq_inv[0]*(-Edp + Vd) + Zdq_inv[1]*(-Eqp + Vq) + Id; */
398:     row[0] = idx+4;
399:     col[0] = idx;         col[1] = idx+1;        col[2] = idx + 2;
400:     val[0] = -Zdq_inv[1]; val[1] = -Zdq_inv[0];  val[2] = Zdq_inv[0]*dVd_ddelta + Zdq_inv[1]*dVq_ddelta;
401:     col[3] = idx + 4; col[4] = net_start+2*gbus[i];                     col[5] = net_start + 2*gbus[i]+1;
402:     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;
403:     MatSetValues(J,1,row,6,col,val,INSERT_VALUES);

405:     /*  fgen[idx+5] = Zdq_inv[2]*(-Edp + Vd) + Zdq_inv[3]*(-Eqp + Vq) + Iq; */
406:     row[0] = idx+5;
407:     col[0] = idx;         col[1] = idx+1;        col[2] = idx + 2;
408:     val[0] = -Zdq_inv[3]; val[1] = -Zdq_inv[2];  val[2] = Zdq_inv[2]*dVd_ddelta + Zdq_inv[3]*dVq_ddelta;
409:     col[3] = idx + 5; col[4] = net_start+2*gbus[i];                     col[5] = net_start + 2*gbus[i]+1;
410:     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;
411:     MatSetValues(J,1,row,6,col,val,INSERT_VALUES);

413:     dIGr_ddelta = Id*PetscCosScalar(delta) - Iq*PetscSinScalar(delta);
414:     dIGi_ddelta = Id*PetscSinScalar(delta) + Iq*PetscCosScalar(delta);
415:     dIGr_dId    = PetscSinScalar(delta);  dIGr_dIq = PetscCosScalar(delta);
416:     dIGi_dId    = -PetscCosScalar(delta); dIGi_dIq = PetscSinScalar(delta);

418:     /* fnet[2*gbus[i]]   -= IGi; */
419:     row[0] = net_start + 2*gbus[i];
420:     col[0] = idx+2;        col[1] = idx + 4;   col[2] = idx + 5;
421:     val[0] = -dIGi_ddelta; val[1] = -dIGi_dId; val[2] = -dIGi_dIq;
422:     MatSetValues(J,1,row,3,col,val,INSERT_VALUES);

424:     /* fnet[2*gbus[i]+1]   -= IGr; */
425:     row[0] = net_start + 2*gbus[i]+1;
426:     col[0] = idx+2;        col[1] = idx + 4;   col[2] = idx + 5;
427:     val[0] = -dIGr_ddelta; val[1] = -dIGr_dId; val[2] = -dIGr_dIq;
428:     MatSetValues(J,1,row,3,col,val,INSERT_VALUES);

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

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

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

437:     row[0] = idx + 6;
438:     col[0] = idx + 6;                     col[1] = idx + 8;
439:     val[0] = (KE[i] + dSE_dEfd)/TE[i];  val[1] = -1/TE[i];
440:     MatSetValues(J,1,row,2,col,val,INSERT_VALUES);

442:     /* Exciter differential equations */

444:     /*    fgen[idx+7] = (RF - KF[i]*Efd/TF[i])/TF[i]; */
445:     row[0] = idx + 7;
446:     col[0] = idx + 6;       col[1] = idx + 7;
447:     val[0] = (-KF[i]/TF[i])/TF[i];  val[1] = 1/TF[i];
448:     MatSetValues(J,1,row,2,col,val,INSERT_VALUES);

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

453:     dVm_dVd = Vd/Vm; dVm_dVq = Vq/Vm;
454:     dVm_dVr = dVm_dVd*dVd_dVr + dVm_dVq*dVq_dVr;
455:     dVm_dVi = dVm_dVd*dVd_dVi + dVm_dVq*dVq_dVi;
456:     row[0]  = idx + 8;
457:     col[0]  = idx + 6;           col[1] = idx + 7; col[2] = idx + 8;
458:     val[0]  = (KA[i]*KF[i]/TF[i])/TA[i]; val[1] = -KA[i]/TA[i];  val[2] = 1/TA[i];
459:     col[3]  = net_start + 2*gbus[i]; col[4] = net_start + 2*gbus[i]+1;
460:     val[3]  = KA[i]*dVm_dVr/TA[i];         val[4] = KA[i]*dVm_dVi/TA[i];
461:     MatSetValues(J,1,row,5,col,val,INSERT_VALUES);
462:     idx     = idx + 9;
463:   }

465:   for (i=0; i<nbus; i++) {
466:     MatGetRow(user->Ybus,2*i,&ncols,&cols,&yvals);
467:     row[0] = net_start + 2*i;
468:     for (k=0; k<ncols; k++) {
469:       col[k] = net_start + cols[k];
470:       val[k] = yvals[k];
471:     }
472:     MatSetValues(J,1,row,ncols,col,val,INSERT_VALUES);
473:     MatRestoreRow(user->Ybus,2*i,&ncols,&cols,&yvals);

475:     MatGetRow(user->Ybus,2*i+1,&ncols,&cols,&yvals);
476:     row[0] = net_start + 2*i+1;
477:     for (k=0; k<ncols; k++) {
478:       col[k] = net_start + cols[k];
479:       val[k] = yvals[k];
480:     }
481:     MatSetValues(J,1,row,ncols,col,val,INSERT_VALUES);
482:     MatRestoreRow(user->Ybus,2*i+1,&ncols,&cols,&yvals);
483:   }

485:   MatAssemblyBegin(J,MAT_FLUSH_ASSEMBLY);
486:   MatAssemblyEnd(J,MAT_FLUSH_ASSEMBLY);

488:   VecGetArray(user->V0,&v0);
489:   for (i=0; i < nload; i++) {
490:     Vr      = xnet[2*lbus[i]]; /* Real part of load bus voltage */
491:     Vi      = xnet[2*lbus[i]+1]; /* Imaginary part of the load bus voltage */
492:     Vm      = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm; Vm4 = Vm2*Vm2;
493:     Vm0     = PetscSqrtScalar(v0[2*lbus[i]]*v0[2*lbus[i]] + v0[2*lbus[i]+1]*v0[2*lbus[i]+1]);
494:     PD      = QD = 0.0;
495:     dPD_dVr = dPD_dVi = dQD_dVr = dQD_dVi = 0.0;
496:     for (k=0; k < ld_nsegsp[i]; k++) {
497:       PD      += ld_alphap[k]*PD0[i]*PetscPowScalar((Vm/Vm0),ld_betap[k]);
498:       dPD_dVr += ld_alphap[k]*ld_betap[k]*PD0[i]*PetscPowScalar((1/Vm0),ld_betap[k])*Vr*PetscPowScalar(Vm,(ld_betap[k]-2));
499:       dPD_dVi += ld_alphap[k]*ld_betap[k]*PD0[i]*PetscPowScalar((1/Vm0),ld_betap[k])*Vi*PetscPowScalar(Vm,(ld_betap[k]-2));
500:     }
501:     for (k=0; k < ld_nsegsq[i]; k++) {
502:       QD      += ld_alphaq[k]*QD0[i]*PetscPowScalar((Vm/Vm0),ld_betaq[k]);
503:       dQD_dVr += ld_alphaq[k]*ld_betaq[k]*QD0[i]*PetscPowScalar((1/Vm0),ld_betaq[k])*Vr*PetscPowScalar(Vm,(ld_betaq[k]-2));
504:       dQD_dVi += ld_alphaq[k]*ld_betaq[k]*QD0[i]*PetscPowScalar((1/Vm0),ld_betaq[k])*Vi*PetscPowScalar(Vm,(ld_betaq[k]-2));
505:     }

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

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

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


517:     /*    fnet[2*lbus[i]]   += IDi; */
518:     row[0] = net_start + 2*lbus[i];
519:     col[0] = net_start + 2*lbus[i];  col[1] = net_start + 2*lbus[i]+1;
520:     val[0] = dIDi_dVr;               val[1] = dIDi_dVi;
521:     MatSetValues(J,1,row,2,col,val,ADD_VALUES);
522:     /*    fnet[2*lbus[i]+1] += IDr; */
523:     row[0] = net_start + 2*lbus[i]+1;
524:     col[0] = net_start + 2*lbus[i];  col[1] = net_start + 2*lbus[i]+1;
525:     val[0] = dIDr_dVr;               val[1] = dIDr_dVi;
526:     MatSetValues(J,1,row,2,col,val,ADD_VALUES);
527:   }
528:   VecRestoreArray(user->V0,&v0);

530:   VecRestoreArray(Xgen,&xgen);
531:   VecRestoreArray(Xnet,&xnet);

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

535:   MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
536:   MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
537:   return(0);
538: }

542: int main(int argc,char **argv)
543: {
544:   EPS            eps;
545:   EPSType        type;
547:   PetscMPIInt    size;
548:   Userctx        user;
549:   PetscViewer    Xview,Ybusview;
550:   Vec            X,Xr,Xi;
551:   Mat            J,Jred=NULL;
552:   IS             is0,is1;
553:   PetscInt       i,*idx2,its,nev,nconv;
554:   PetscReal      error,re,im;
555:   PetscScalar    kr,ki;
556:   PetscBool      terse;

558:   SlepcInitialize(&argc,&argv,NULL,help);
559:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
560:   if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs");
561:   /* show detailed info unless -terse option is given by user */
562:   PetscOptionsHasName(NULL,NULL,"-terse",&terse);

564:   user.neqs_gen   = 9*ngen; /* # eqs. for generator subsystem */
565:   user.neqs_net   = 2*nbus; /* # eqs. for network subsystem   */
566:   user.neqs_pgrid = user.neqs_gen + user.neqs_net;
567:   PetscPrintf(PETSC_COMM_WORLD,"\nStability analysis in a network with %D buses and %D generators\n\n",nbus,ngen);

569:   /* Create indices for differential and algebraic equations */
570:   PetscMalloc1(7*ngen,&idx2);
571:   for (i=0; i<ngen; i++) {
572:     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;
573:     idx2[7*i+4] = 9*i+6; idx2[7*i+5] = 9*i+7; idx2[7*i+6] = 9*i+8;
574:   }
575:   ISCreateGeneral(PETSC_COMM_WORLD,7*ngen,idx2,PETSC_COPY_VALUES,&user.is_diff);
576:   ISComplement(user.is_diff,0,user.neqs_pgrid,&user.is_alg);
577:   PetscFree(idx2);

579:   /* Read initial voltage vector and Ybus */
580:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"X.bin",FILE_MODE_READ,&Xview);
581:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"Ybus.bin",FILE_MODE_READ,&Ybusview);

583:   VecCreate(PETSC_COMM_WORLD,&user.V0);
584:   VecSetSizes(user.V0,PETSC_DECIDE,user.neqs_net);
585:   VecLoad(user.V0,Xview);

587:   MatCreate(PETSC_COMM_WORLD,&user.Ybus);
588:   MatSetSizes(user.Ybus,PETSC_DECIDE,PETSC_DECIDE,user.neqs_net,user.neqs_net);
589:   MatSetType(user.Ybus,MATBAIJ);
590:   /*  MatSetBlockSize(user.Ybus,2); */
591:   MatLoad(user.Ybus,Ybusview);

593:   PetscViewerDestroy(&Xview);
594:   PetscViewerDestroy(&Ybusview);

596:   /* Create DMs for generator and network subsystems */
597:   DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,user.neqs_gen,1,1,NULL,&user.dmgen);
598:   DMSetOptionsPrefix(user.dmgen,"dmgen_");
599:   DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,user.neqs_net,1,1,NULL,&user.dmnet);
600:   DMSetOptionsPrefix(user.dmnet,"dmnet_");
601:   /* Create a composite DM packer and add the two DMs */
602:   DMCompositeCreate(PETSC_COMM_WORLD,&user.dmpgrid);
603:   DMSetOptionsPrefix(user.dmpgrid,"pgrid_");
604:   DMCompositeAddDM(user.dmpgrid,user.dmgen);
605:   DMCompositeAddDM(user.dmpgrid,user.dmnet);

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

609:   MatCreate(PETSC_COMM_WORLD,&J);
610:   MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,user.neqs_pgrid,user.neqs_pgrid);
611:   MatSetFromOptions(J);
612:   PreallocateJacobian(J,&user);

614:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
615:      Set initial conditions
616:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
617:   SetInitialGuess(X,&user);

619:   /* Form Jacobian */
620:   ResidualJacobian(X,J,(void*)&user);
621:   MatScale(J,-1);
622:   is0 = user.is_diff;
623:   is1 = user.is_alg;

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

627:   if (!terse) {
628:     MatView(Jred,NULL);
629:   }

631:   MatCreateVecs(Jred,NULL,&Xr);
632:   MatCreateVecs(Jred,NULL,&Xi);

634:   /* Create the eigensolver and set the various options */
635:   EPSCreate(PETSC_COMM_WORLD,&eps);
636:   EPSSetOperators(eps,Jred,NULL);
637:   EPSSetProblemType(eps,EPS_NHEP);
638:   EPSSetFromOptions(eps);
639:   
640:   /* Solve the eigenvalue problem */
641:   EPSSolve(eps);

643:   EPSGetIterationNumber(eps,&its);
644:   PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the eigensolver: %D\n",its);
645:   EPSGetType(eps,&type);
646:   PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n", type);
647:   EPSGetDimensions(eps,&nev,NULL,NULL);
648:   PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);

650:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
651:                     Display solution and clean up
652:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
653:   if (terse) {
654:     EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
655:   } else {
656:     /* Get number of converged approximate eigenpairs */
657:     EPSGetConverged(eps,&nconv);
658:     PetscPrintf(PETSC_COMM_WORLD," Number of converged eigenpairs: %D\n\n",nconv);
659:   
660:     if (nconv>0) {
661:       /* Display eigenvalues and relative errors */
662:       PetscPrintf(PETSC_COMM_WORLD,
663:            "           k          ||Ax-kx||/||kx||\n"
664:            "   ----------------- ------------------\n");
665:   
666:       for (i=0;i<nconv;i++) {
667:         /* Get converged eigenpairs: i-th eigenvalue is stored in kr (real part) and
668:           ki (imaginary part) */
669:         EPSGetEigenpair(eps,i,&kr,&ki,Xr,Xi);
670:         /* Compute the relative error associated to each eigenpair */
671:         EPSComputeError(eps,i,EPS_ERROR_RELATIVE,&error);

673: #if defined(PETSC_USE_COMPLEX)
674:         re = PetscRealPart(kr);
675:         im = PetscImaginaryPart(kr);
676: #else
677:         re = kr;
678:         im = ki;
679: #endif
680:         if (im!=0.0) {
681:           PetscPrintf(PETSC_COMM_WORLD," %9f%+9fi %12g\n",(double)re,(double)im,(double)error);
682:         } else {
683:           PetscPrintf(PETSC_COMM_WORLD,"   %12f       %12g\n",(double)re,(double)error);
684:         }
685:       }
686:       PetscPrintf(PETSC_COMM_WORLD,"\n");
687:     }
688:   }

690:   /* Free work space */
691:   EPSDestroy(&eps);
692:   MatDestroy(&J);
693:   MatDestroy(&Jred);
694:   MatDestroy(&user.Ybus);
695:   VecDestroy(&X);
696:   VecDestroy(&Xr);
697:   VecDestroy(&Xi);
698:   VecDestroy(&user.V0);
699:   DMDestroy(&user.dmgen);
700:   DMDestroy(&user.dmnet);
701:   DMDestroy(&user.dmpgrid);
702:   ISDestroy(&user.is_diff);
703:   ISDestroy(&user.is_alg);
704:   SlepcFinalize();
705:   return ierr;
706: }