Actual source code: ex9bus.c

petsc-3.8.3 2017-12-09
Report Typos and Errors

  2: static char help[] = "Power grid stability analysis of WECC 9 bus system.\n\
  3: This example is based on the 9-bus (node) example given in the book Power\n\
  4: Systems Dynamics and Stability (Chapter 7) by P. Sauer and M. A. Pai.\n\
  5: The power grid in this example consists of 9 buses (nodes), 3 generators,\n\
  6: 3 loads, and 9 transmission lines. The network equations are written\n\
  7: in current balance form using rectangular coordiantes.\n\n";

  9: /*
 10:    The equations for the stability analysis are described by the DAE

 12:    \dot{x} = f(x,y,t)
 13:      0     = g(x,y,t)

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

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

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

 26:    where:
 27:     I(x,y) is the current injected from generators and loads.
 28:       Y    is the admittance matrix, and
 29:       V    is the voltage vector
 30: */

 32: /*
 33:    Include "petscts.h" so that we can use TS solvers.  Note that this
 34:    file automatically includes:
 35:      petscsys.h       - base PETSc routines   petscvec.h - vectors
 36:      petscmat.h - matrices
 37:      petscis.h     - index sets            petscksp.h - Krylov subspace methods
 38:      petscviewer.h - viewers               petscpc.h  - preconditioners
 39:      petscksp.h   - linear solvers
 40: */
 41: #include <petscts.h>
 42: #include <petscdm.h>
 43: #include <petscdmda.h>
 44: #include <petscdmcomposite.h>

 46: #define freq 60
 47: #define w_s (2*PETSC_PI*freq)

 49: /* Sizes and indices */
 50: const PetscInt nbus    = 9; /* Number of network buses */
 51: const PetscInt ngen    = 3; /* Number of generators */
 52: const PetscInt nload   = 3; /* Number of loads */
 53: const PetscInt gbus[3] = {0,1,2}; /* Buses at which generators are incident */
 54: const PetscInt lbus[3] = {4,5,7}; /* Buses at which loads are incident */

 56: /* Generator real and reactive powers (found via loadflow) */
 57: const PetscScalar PG[3] = {0.716786142395021,1.630000000000000,0.850000000000000};
 58: const PetscScalar QG[3] = {0.270702180178785,0.066120127797275,-0.108402221791588};
 59: /* Generator constants */
 60: const PetscScalar H[3]    = {23.64,6.4,3.01};   /* Inertia constant */
 61: const PetscScalar Rs[3]   = {0.0,0.0,0.0}; /* Stator Resistance */
 62: const PetscScalar Xd[3]   = {0.146,0.8958,1.3125};  /* d-axis reactance */
 63: const PetscScalar Xdp[3]  = {0.0608,0.1198,0.1813}; /* d-axis transient reactance */
 64: const PetscScalar Xq[3]   = {0.4360,0.8645,1.2578}; /* q-axis reactance Xq(1) set to 0.4360, value given in text 0.0969 */
 65: const PetscScalar Xqp[3]  = {0.0969,0.1969,0.25}; /* q-axis transient reactance */
 66: const PetscScalar Td0p[3] = {8.96,6.0,5.89}; /* d-axis open circuit time constant */
 67: const PetscScalar Tq0p[3] = {0.31,0.535,0.6}; /* q-axis open circuit time constant */
 68: PetscScalar M[3]; /* M = 2*H/w_s */
 69: PetscScalar D[3]; /* D = 0.1*M */

 71: PetscScalar TM[3]; /* Mechanical Torque */
 72: /* Exciter system constants */
 73: const PetscScalar KA[3] = {20.0,20.0,20.0};  /* Voltage regulartor gain constant */
 74: const PetscScalar TA[3] = {0.2,0.2,0.2};     /* Voltage regulator time constant */
 75: const PetscScalar KE[3] = {1.0,1.0,1.0};     /* Exciter gain constant */
 76: const PetscScalar TE[3] = {0.314,0.314,0.314}; /* Exciter time constant */
 77: const PetscScalar KF[3] = {0.063,0.063,0.063};  /* Feedback stabilizer gain constant */
 78: const PetscScalar TF[3] = {0.35,0.35,0.35};    /* Feedback stabilizer time constant */
 79: const PetscScalar k1[3] = {0.0039,0.0039,0.0039};
 80: const PetscScalar k2[3] = {1.555,1.555,1.555};  /* k1 and k2 for calculating the saturation function SE = k1*exp(k2*Efd) */
 81: const PetscScalar VRMIN[3] = {-4.0,-4.0,-4.0};
 82: const PetscScalar VRMAX[3] = {7.0,7.0,7.0};
 83: PetscInt VRatmin[3];
 84: PetscInt VRatmax[3];

 86: PetscScalar Vref[3];
 87: /* Load constants
 88:   We use a composite load model that describes the load and reactive powers at each time instant as follows
 89:   P(t) = \sum\limits_{i=0}^ld_nsegsp \ld_alphap_i*P_D0(\frac{V_m(t)}{V_m0})^\ld_betap_i
 90:   Q(t) = \sum\limits_{i=0}^ld_nsegsq \ld_alphaq_i*Q_D0(\frac{V_m(t)}{V_m0})^\ld_betaq_i
 91:   where
 92:     ld_nsegsp,ld_nsegsq - Number of individual load models for real and reactive power loads
 93:     ld_alphap,ld_alphap - Percentage contribution (weights) or loads
 94:     P_D0                - Real power load
 95:     Q_D0                - Reactive power load
 96:     V_m(t)              - Voltage magnitude at time t
 97:     V_m0                - Voltage magnitude at t = 0
 98:     ld_betap, ld_betaq  - exponents describing the load model for real and reactive part

100:     Note: All loads have the same characteristic currently.
101: */
102: const PetscScalar PD0[3] = {1.25,0.9,1.0};
103: const PetscScalar QD0[3] = {0.5,0.3,0.35};
104: const PetscInt    ld_nsegsp[3] = {3,3,3};
105: const PetscScalar ld_alphap[3] = {1.0,0.0,0.0};
106: const PetscScalar ld_betap[3]  = {2.0,1.0,0.0};
107: const PetscInt    ld_nsegsq[3] = {3,3,3};
108: const PetscScalar ld_alphaq[3] = {1.0,0.0,0.0};
109: const PetscScalar ld_betaq[3]  = {2.0,1.0,0.0};

111: typedef struct {
112:   DM          dmgen, dmnet; /* DMs to manage generator and network subsystem */
113:   DM          dmpgrid; /* Composite DM to manage the entire power grid */
114:   Mat         Ybus; /* Network admittance matrix */
115:   Vec         V0;  /* Initial voltage vector (Power flow solution) */
116:   PetscReal   tfaulton,tfaultoff; /* Fault on and off times */
117:   PetscInt    faultbus; /* Fault bus */
118:   PetscScalar Rfault;
119:   PetscReal   t0,tmax;
120:   PetscInt    neqs_gen,neqs_net,neqs_pgrid;
121:   Mat         Sol; /* Matrix to save solution at each time step */
122:   PetscInt    stepnum;
123:   PetscReal   t;
124:   SNES        snes_alg;
125:   IS          is_diff; /* indices for differential equations */
126:   IS          is_alg; /* indices for algebraic equations */
127:   PetscBool   setisdiff; /* TS computes truncation error based only on the differential variables */
128:   PetscBool   semiexplicit; /* If the flag is set then a semi-explicit method is used using TSRK */
129: } Userctx;

131: /*
132:   The first two events are for fault on and off, respectively. The following events are
133:   to check the min/max limits on the state variable VR. A non windup limiter is used for
134:   the VR limits.
135: */
136: PetscErrorCode EventFunction(TS ts,PetscReal t,Vec X,PetscScalar *fvalue,void *ctx)
137: {
138:   Userctx        *user=(Userctx*)ctx;
139:   Vec            Xgen,Xnet;
140:   PetscInt       i,idx=0;
141:   const PetscScalar *xgen,*xnet;
143:   PetscScalar    Efd,RF,VR,Vr,Vi,Vm;


147:   DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
148:   DMCompositeScatter(user->dmpgrid,X,Xgen,Xnet);

150:   VecGetArrayRead(Xgen,&xgen);
151:   VecGetArrayRead(Xnet,&xnet);

153:   /* Event for fault-on time */
154:   fvalue[0] = t - user->tfaulton;
155:   /* Event for fault-off time */
156:   fvalue[1] = t - user->tfaultoff;

158:   for (i=0; i < ngen; i++) {
159:     Efd   = xgen[idx+6];
160:     RF    = xgen[idx+7];
161:     VR    = xgen[idx+8];

163:     Vr = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */
164:     Vi = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */
165:     Vm = PetscSqrtScalar(Vr*Vr + Vi*Vi);

167:     if (!VRatmax[i]) {
168:       fvalue[2+2*i] = VRMAX[i] - VR;
169:     } else {
170:       fvalue[2+2*i] = (VR - KA[i]*RF + KA[i]*KF[i]*Efd/TF[i] - KA[i]*(Vref[i] - Vm))/TA[i];
171:     }
172:     if (!VRatmin[i]) {
173:       fvalue[2+2*i+1] = VRMIN[i] - VR;
174:     } else {
175:       fvalue[2+2*i+1] = (VR - KA[i]*RF + KA[i]*KF[i]*Efd/TF[i] - KA[i]*(Vref[i] - Vm))/TA[i];
176:     }
177:     idx = idx+9;
178:   }
179:   VecRestoreArrayRead(Xgen,&xgen);
180:   VecRestoreArrayRead(Xnet,&xnet);

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

184:   return(0);
185: }

187: PetscErrorCode PostEventFunction(TS ts,PetscInt nevents,PetscInt event_list[],PetscReal t,Vec X,PetscBool forwardsolve,void* ctx)
188: {
189:   Userctx *user=(Userctx*)ctx;
190:   Vec      Xgen,Xnet;
191:   PetscScalar *xgen,*xnet;
192:   PetscInt row_loc,col_loc;
193:   PetscScalar val;
195:   PetscInt i,idx=0,event_num;
196:   PetscScalar fvalue;
197:   PetscScalar Efd, RF, VR;
198:   PetscScalar Vr,Vi,Vm;
199: 

202:   DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
203:   DMCompositeScatter(user->dmpgrid,X,Xgen,Xnet);

205:   VecGetArray(Xgen,&xgen);
206:   VecGetArray(Xnet,&xnet);

208:   for (i=0; i < nevents; i++) {
209:     if (event_list[i] == 0) {
210:       /* Apply disturbance - resistive fault at user->faultbus */
211:       /* This is done by adding shunt conductance to the diagonal location
212:          in the Ybus matrix */
213:       row_loc = 2*user->faultbus; col_loc = 2*user->faultbus+1; /* Location for G */
214:       val     = 1/user->Rfault;
215:       MatSetValues(user->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);
216:       row_loc = 2*user->faultbus+1; col_loc = 2*user->faultbus; /* Location for G */
217:       val     = 1/user->Rfault;
218:       MatSetValues(user->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);
219: 
220:       MatAssemblyBegin(user->Ybus,MAT_FINAL_ASSEMBLY);
221:       MatAssemblyEnd(user->Ybus,MAT_FINAL_ASSEMBLY);
222: 
223:       /* Solve the algebraic equations */
224:       SNESSolve(user->snes_alg,NULL,X);
225:     } else if(event_list[i] == 1) {
226:       /* Remove the fault */
227:       row_loc = 2*user->faultbus; col_loc = 2*user->faultbus+1;
228:       val     = -1/user->Rfault;
229:       MatSetValues(user->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);
230:       row_loc = 2*user->faultbus+1; col_loc = 2*user->faultbus;
231:       val     = -1/user->Rfault;
232:       MatSetValues(user->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);
233: 
234:       MatAssemblyBegin(user->Ybus,MAT_FINAL_ASSEMBLY);
235:       MatAssemblyEnd(user->Ybus,MAT_FINAL_ASSEMBLY);
236: 
237:       /* Solve the algebraic equations */
238:       SNESSolve(user->snes_alg,NULL,X);

240:       /* Check the VR derivatives and reset flags if needed */
241:       for (i=0; i < ngen; i++) {
242:         Efd   = xgen[idx+6];
243:         RF    = xgen[idx+7];
244:         VR    = xgen[idx+8];

246:         Vr = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */
247:         Vi = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */
248:         Vm = PetscSqrtScalar(Vr*Vr + Vi*Vi);

250:         if (VRatmax[i]) {
251:           fvalue = (VR - KA[i]*RF + KA[i]*KF[i]*Efd/TF[i] - KA[i]*(Vref[i] - Vm))/TA[i];
252:           if (fvalue < 0) {
253:             VRatmax[i] = 0;
254:             PetscPrintf(PETSC_COMM_SELF,"VR[%d]: dVR_dt went negative on fault clearing at time %g\n",i,t);
255:           }
256:         }
257:         if (VRatmin[i]) {
258:           fvalue =  (VR - KA[i]*RF + KA[i]*KF[i]*Efd/TF[i] - KA[i]*(Vref[i] - Vm))/TA[i];

260:           if(fvalue > 0) {
261:             VRatmin[i] = 0;
262:             PetscPrintf(PETSC_COMM_SELF,"VR[%d]: dVR_dt went positive on fault clearing at time %g\n",i,t);
263:           }
264:         }
265:         idx = idx+9;
266:       }
267:     } else {
268:       idx = (event_list[i]-2)/2;
269:       event_num = (event_list[i]-2)%2;
270:       if (event_num == 0) { /* Max VR */
271:         if (!VRatmax[idx]) {
272:           VRatmax[idx] = 1;
273:           PetscPrintf(PETSC_COMM_SELF,"VR[%d]: hit upper limit at time %g\n",idx,t);
274:         }
275:         else {
276:           VRatmax[idx] = 0;
277:           PetscPrintf(PETSC_COMM_SELF,"VR[%d]: freeing variable as dVR_dt is negative at time %g\n",idx,t);
278:         }
279:       } else {
280:         if (!VRatmin[idx]) {
281:           VRatmin[idx] = 1;
282:           PetscPrintf(PETSC_COMM_SELF,"VR[%d]: hit lower limit at time %g\n",idx,t);
283:         }
284:         else {
285:           VRatmin[idx] = 0;
286:           PetscPrintf(PETSC_COMM_SELF,"VR[%d]: freeing variable as dVR_dt is positive at time %g\n",idx,t);
287:         }
288:       }
289:     }
290:   }

292:   VecRestoreArray(Xgen,&xgen);
293:   VecRestoreArray(Xnet,&xnet);

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

297:   return(0);
298: }

300: /* Converts from machine frame (dq) to network (phase a real,imag) reference frame */
301: PetscErrorCode dq2ri(PetscScalar Fd,PetscScalar Fq,PetscScalar delta,PetscScalar *Fr, PetscScalar *Fi)
302: {
304:   *Fr =  Fd*PetscSinScalar(delta) + Fq*PetscCosScalar(delta);
305:   *Fi = -Fd*PetscCosScalar(delta) + Fq*PetscSinScalar(delta);
306:   return(0);
307: }

309: /* Converts from network frame ([phase a real,imag) to machine (dq) reference frame */
310: PetscErrorCode ri2dq(PetscScalar Fr,PetscScalar Fi,PetscScalar delta,PetscScalar *Fd, PetscScalar *Fq)
311: {
313:   *Fd =  Fr*PetscSinScalar(delta) - Fi*PetscCosScalar(delta);
314:   *Fq =  Fr*PetscCosScalar(delta) + Fi*PetscSinScalar(delta);
315:   return(0);
316: }

318: /* Saves the solution at each time to a matrix */
319: PetscErrorCode SaveSolution(TS ts)
320: {
321:   PetscErrorCode    ierr;
322:   Userctx           *user;
323:   Vec               X;
324:   const PetscScalar *x;
325:   PetscScalar       *mat;
326:   PetscInt          idx;
327:   PetscReal         t;

330:   TSGetApplicationContext(ts,&user);
331:   TSGetTime(ts,&t);
332:   TSGetSolution(ts,&X);
333:   idx      = user->stepnum*(user->neqs_pgrid+1);
334:   MatDenseGetArray(user->Sol,&mat);
335:   VecGetArrayRead(X,&x);
336:   mat[idx] = t;
337:   PetscMemcpy(mat+idx+1,x,user->neqs_pgrid*sizeof(PetscScalar));
338:   MatDenseRestoreArray(user->Sol,&mat);
339:   VecRestoreArrayRead(X,&x);
340:   user->stepnum++;
341:   return(0);
342: }

344: PetscErrorCode SetInitialGuess(Vec X,Userctx *user)
345: {
347:   Vec            Xgen,Xnet;
348:   PetscScalar    *xgen,*xnet;
349:   PetscInt       i,idx=0;
350:   PetscScalar    Vr,Vi,IGr,IGi,Vm,Vm2;
351:   PetscScalar    Eqp,Edp,delta;
352:   PetscScalar    Efd,RF,VR; /* Exciter variables */
353:   PetscScalar    Id,Iq;  /* Generator dq axis currents */
354:   PetscScalar    theta,Vd,Vq,SE;

357:   M[0] = 2*H[0]/w_s; M[1] = 2*H[1]/w_s; M[2] = 2*H[2]/w_s;
358:   D[0] = 0.1*M[0]; D[1] = 0.1*M[1]; D[2] = 0.1*M[2];

360:   DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);

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

365:   /* Generator subsystem initialization */
366:   VecGetArray(Xgen,&xgen);
367:   VecGetArray(Xnet,&xnet);

369:   for (i=0; i < ngen; i++) {
370:     Vr  = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */
371:     Vi  = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */
372:     Vm  = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm;
373:     IGr = (Vr*PG[i] + Vi*QG[i])/Vm2;
374:     IGi = (Vi*PG[i] - Vr*QG[i])/Vm2;

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

378:     theta = PETSC_PI/2.0 - delta;

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

383:     Vd = Vr*PetscCosScalar(theta) - Vi*PetscSinScalar(theta);
384:     Vq = Vr*PetscSinScalar(theta) + Vi*PetscCosScalar(theta);

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

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

391:     /* The generator variables are ordered as [Eqp,Edp,delta,w,Id,Iq] */
392:     xgen[idx]   = Eqp;
393:     xgen[idx+1] = Edp;
394:     xgen[idx+2] = delta;
395:     xgen[idx+3] = w_s;

397:     idx = idx + 4;

399:     xgen[idx]   = Id;
400:     xgen[idx+1] = Iq;

402:     idx = idx + 2;

404:     /* Exciter */
405:     Efd = Eqp + (Xd[i] - Xdp[i])*Id;
406:     SE  = k1[i]*PetscExpScalar(k2[i]*Efd);
407:     VR  =  KE[i]*Efd + SE;
408:     RF  =  KF[i]*Efd/TF[i];

410:     xgen[idx]   = Efd;
411:     xgen[idx+1] = RF;
412:     xgen[idx+2] = VR;

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

416:     VRatmin[i] = VRatmax[i] = 0;

418:     idx = idx + 3;
419:   }

421:   VecRestoreArray(Xgen,&xgen);
422:   VecRestoreArray(Xnet,&xnet);

424:   /* VecView(Xgen,0); */
425:   DMCompositeGather(user->dmpgrid,INSERT_VALUES,X,Xgen,Xnet);
426:   DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);
427:   return(0);
428: }

430: /* Computes F = [f(x,y);g(x,y)] */
431: PetscErrorCode ResidualFunction(Vec X, Vec F, Userctx *user)
432: {
434:   Vec            Xgen,Xnet,Fgen,Fnet;
435:   PetscScalar    *xgen,*xnet,*fgen,*fnet;
436:   PetscInt       i,idx=0;
437:   PetscScalar    Vr,Vi,Vm,Vm2;
438:   PetscScalar    Eqp,Edp,delta,w; /* Generator variables */
439:   PetscScalar    Efd,RF,VR; /* Exciter variables */
440:   PetscScalar    Id,Iq;  /* Generator dq axis currents */
441:   PetscScalar    Vd,Vq,SE;
442:   PetscScalar    IGr,IGi,IDr,IDi;
443:   PetscScalar    Zdq_inv[4],det;
444:   PetscScalar    PD,QD,Vm0,*v0;
445:   PetscInt       k;

448:   VecZeroEntries(F);
449:   DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
450:   DMCompositeGetLocalVectors(user->dmpgrid,&Fgen,&Fnet);
451:   DMCompositeScatter(user->dmpgrid,X,Xgen,Xnet);
452:   DMCompositeScatter(user->dmpgrid,F,Fgen,Fnet);

454:   /* Network current balance residual IG + Y*V + IL = 0. Only YV is added here.
455:      The generator current injection, IG, and load current injection, ID are added later
456:   */
457:   /* Note that the values in Ybus are stored assuming the imaginary current balance
458:      equation is ordered first followed by real current balance equation for each bus.
459:      Thus imaginary current contribution goes in location 2*i, and
460:      real current contribution in 2*i+1
461:   */
462:   MatMult(user->Ybus,Xnet,Fnet);

464:   VecGetArray(Xgen,&xgen);
465:   VecGetArray(Xnet,&xnet);
466:   VecGetArray(Fgen,&fgen);
467:   VecGetArray(Fnet,&fnet);

469:   /* Generator subsystem */
470:   for (i=0; i < ngen; i++) {
471:     Eqp   = xgen[idx];
472:     Edp   = xgen[idx+1];
473:     delta = xgen[idx+2];
474:     w     = xgen[idx+3];
475:     Id    = xgen[idx+4];
476:     Iq    = xgen[idx+5];
477:     Efd   = xgen[idx+6];
478:     RF    = xgen[idx+7];
479:     VR    = xgen[idx+8];

481:     /* Generator differential equations */
482:     fgen[idx]   = (-Eqp - (Xd[i] - Xdp[i])*Id + Efd)/Td0p[i];
483:     fgen[idx+1] = (-Edp + (Xq[i] - Xqp[i])*Iq)/Tq0p[i];
484:     fgen[idx+2] = w - w_s;
485:     fgen[idx+3] = (TM[i] - Edp*Id - Eqp*Iq - (Xqp[i] - Xdp[i])*Id*Iq - D[i]*(w - w_s))/M[i];

487:     Vr = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */
488:     Vi = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */

490:     ri2dq(Vr,Vi,delta,&Vd,&Vq);
491:     /* Algebraic equations for stator currents */
492:     det = Rs[i]*Rs[i] + Xdp[i]*Xqp[i];

494:     Zdq_inv[0] = Rs[i]/det;
495:     Zdq_inv[1] = Xqp[i]/det;
496:     Zdq_inv[2] = -Xdp[i]/det;
497:     Zdq_inv[3] = Rs[i]/det;

499:     fgen[idx+4] = Zdq_inv[0]*(-Edp + Vd) + Zdq_inv[1]*(-Eqp + Vq) + Id;
500:     fgen[idx+5] = Zdq_inv[2]*(-Edp + Vd) + Zdq_inv[3]*(-Eqp + Vq) + Iq;

502:     /* Add generator current injection to network */
503:     dq2ri(Id,Iq,delta,&IGr,&IGi);

505:     fnet[2*gbus[i]]   -= IGi;
506:     fnet[2*gbus[i]+1] -= IGr;

508:     Vm = PetscSqrtScalar(Vd*Vd + Vq*Vq);

510:     SE = k1[i]*PetscExpScalar(k2[i]*Efd);

512:     /* Exciter differential equations */
513:     fgen[idx+6] = (-KE[i]*Efd - SE + VR)/TE[i];
514:     fgen[idx+7] = (-RF + KF[i]*Efd/TF[i])/TF[i];
515:     if(VRatmax[i]) fgen[idx+8] = VR - VRMAX[i];
516:     else if(VRatmin[i]) fgen[idx+8] = VRMIN[i] - VR;
517:     else fgen[idx+8] = (-VR + KA[i]*RF - KA[i]*KF[i]*Efd/TF[i] + KA[i]*(Vref[i] - Vm))/TA[i];

519:     idx = idx + 9;
520:   }

522:   VecGetArray(user->V0,&v0);
523:   for (i=0; i < nload; i++) {
524:     Vr  = xnet[2*lbus[i]]; /* Real part of load bus voltage */
525:     Vi  = xnet[2*lbus[i]+1]; /* Imaginary part of the load bus voltage */
526:     Vm  = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm;
527:     Vm0 = PetscSqrtScalar(v0[2*lbus[i]]*v0[2*lbus[i]] + v0[2*lbus[i]+1]*v0[2*lbus[i]+1]);
528:     PD  = QD = 0.0;
529:     for (k=0; k < ld_nsegsp[i]; k++) PD += ld_alphap[k]*PD0[i]*PetscPowScalar((Vm/Vm0),ld_betap[k]);
530:     for (k=0; k < ld_nsegsq[i]; k++) QD += ld_alphaq[k]*QD0[i]*PetscPowScalar((Vm/Vm0),ld_betaq[k]);

532:     /* Load currents */
533:     IDr = (PD*Vr + QD*Vi)/Vm2;
534:     IDi = (-QD*Vr + PD*Vi)/Vm2;

536:     fnet[2*lbus[i]]   += IDi;
537:     fnet[2*lbus[i]+1] += IDr;
538:   }
539:   VecRestoreArray(user->V0,&v0);

541:   VecRestoreArray(Xgen,&xgen);
542:   VecRestoreArray(Xnet,&xnet);
543:   VecRestoreArray(Fgen,&fgen);
544:   VecRestoreArray(Fnet,&fnet);

546:   DMCompositeGather(user->dmpgrid,INSERT_VALUES,F,Fgen,Fnet);
547:   DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);
548:   DMCompositeRestoreLocalVectors(user->dmpgrid,&Fgen,&Fnet);
549:   return(0);
550: }

552: /*   f(x,y)
553:      g(x,y)
554:  */
555: PetscErrorCode RHSFunction(TS ts,PetscReal t, Vec X, Vec F, void *ctx)
556: {
558:   Userctx        *user=(Userctx*)ctx;

561:   user->t = t;
562:   ResidualFunction(X,F,user);
563:   return(0);
564: }

566: /* f(x,y) - \dot{x}
567:      g(x,y) = 0
568:  */
569: PetscErrorCode IFunction(TS ts,PetscReal t, Vec X, Vec Xdot, Vec F, void *ctx)
570: {
572:   PetscScalar    *f,*xdot;
573:   PetscInt       i;


577:   RHSFunction(ts,t,X,F,ctx);
578:   VecGetArray(F,&f);
579:   VecGetArray(Xdot,&xdot);
580:   for (i=0;i < ngen;i++) {
581:     f[9*i]   -= xdot[9*i];
582:     f[9*i+1] -= xdot[9*i+1];
583:     f[9*i+2] -= xdot[9*i+2];
584:     f[9*i+3] -= xdot[9*i+3];
585:     f[9*i+6] -= xdot[9*i+6];
586:     f[9*i+7] -= xdot[9*i+7];
587:     f[9*i+8] -= xdot[9*i+8];
588:   }
589:   VecRestoreArray(F,&f);
590:   VecRestoreArray(Xdot,&xdot);
591:   return(0);
592: }

594: /* This function is used for solving the algebraic system only during fault on and
595:    off times. It computes the entire F and then zeros out the part corresponding to
596:    differential equations
597:  F = [0;g(y)];
598: */
599: PetscErrorCode AlgFunction(SNES snes, Vec X, Vec F, void *ctx)
600: {
602:   Userctx        *user=(Userctx*)ctx;
603:   PetscScalar    *f;
604:   PetscInt       i;

607:   ResidualFunction(X,F,user);
608:   VecGetArray(F,&f);
609:   for (i=0; i < ngen; i++) {
610:     f[9*i]   = 0;
611:     f[9*i+1] = 0;
612:     f[9*i+2] = 0;
613:     f[9*i+3] = 0;
614:     f[9*i+6] = 0;
615:     f[9*i+7] = 0;
616:     f[9*i+8] = 0;
617:   }
618:   VecRestoreArray(F,&f);
619:   return(0);
620: }

622: PetscErrorCode PostStage(TS ts, PetscReal t, PetscInt i, Vec *X)
623: {
625:   Userctx        *user;

628:   TSGetApplicationContext(ts,&user);
629:   SNESSolve(user->snes_alg,NULL,X[i]);
630:   return(0);
631: }

633: PetscErrorCode PostEvaluate(TS ts)
634: {
636:   Userctx        *user;
637:   Vec            X;

640:   TSGetApplicationContext(ts,&user);
641:   TSGetSolution(ts,&X);
642:   SNESSolve(user->snes_alg,NULL,X);
643:   return(0);
644: }


647: PetscErrorCode PreallocateJacobian(Mat J, Userctx *user)
648: {
650:   PetscInt       *d_nnz;
651:   PetscInt       i,idx=0,start=0;
652:   PetscInt       ncols;

655:   PetscMalloc1(user->neqs_pgrid,&d_nnz);
656:   for (i=0; i<user->neqs_pgrid; i++) d_nnz[i] = 0;
657:   /* Generator subsystem */
658:   for (i=0; i < ngen; i++) {

660:     d_nnz[idx]   += 3;
661:     d_nnz[idx+1] += 2;
662:     d_nnz[idx+2] += 2;
663:     d_nnz[idx+3] += 5;
664:     d_nnz[idx+4] += 6;
665:     d_nnz[idx+5] += 6;

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

670:     d_nnz[idx+6] += 2;
671:     d_nnz[idx+7] += 2;
672:     d_nnz[idx+8] += 5;

674:     idx = idx + 9;
675:   }

677:   start = user->neqs_gen;

679:   for (i=0; i < nbus; i++) {
680:     MatGetRow(user->Ybus,2*i,&ncols,NULL,NULL);
681:     d_nnz[start+2*i]   += ncols;
682:     d_nnz[start+2*i+1] += ncols;
683:     MatRestoreRow(user->Ybus,2*i,&ncols,NULL,NULL);
684:   }

686:   MatSeqAIJSetPreallocation(J,0,d_nnz);

688:   PetscFree(d_nnz);
689:   return(0);
690: }

692: /*
693:    J = [df_dx, df_dy
694:         dg_dx, dg_dy]
695: */
696: PetscErrorCode ResidualJacobian(Vec X,Mat J,Mat B,void *ctx)
697: {
699:   Userctx        *user=(Userctx*)ctx;
700:   Vec            Xgen,Xnet;
701:   PetscScalar    *xgen,*xnet;
702:   PetscInt       i,idx=0;
703:   PetscScalar    Vr,Vi,Vm,Vm2;
704:   PetscScalar    Eqp,Edp,delta; /* Generator variables */
705:   PetscScalar    Efd;
706:   PetscScalar    Id,Iq;  /* Generator dq axis currents */
707:   PetscScalar    Vd,Vq;
708:   PetscScalar    val[10];
709:   PetscInt       row[2],col[10];
710:   PetscInt       net_start=user->neqs_gen;
711:   PetscScalar    Zdq_inv[4],det;
712:   PetscScalar    dVd_dVr,dVd_dVi,dVq_dVr,dVq_dVi,dVd_ddelta,dVq_ddelta;
713:   PetscScalar    dIGr_ddelta,dIGi_ddelta,dIGr_dId,dIGr_dIq,dIGi_dId,dIGi_dIq;
714:   PetscScalar    dSE_dEfd;
715:   PetscScalar    dVm_dVd,dVm_dVq,dVm_dVr,dVm_dVi;
716:   PetscInt          ncols;
717:   const PetscInt    *cols;
718:   const PetscScalar *yvals;
719:   PetscInt          k;
720:   PetscScalar PD,QD,Vm0,*v0,Vm4;
721:   PetscScalar dPD_dVr,dPD_dVi,dQD_dVr,dQD_dVi;
722:   PetscScalar dIDr_dVr,dIDr_dVi,dIDi_dVr,dIDi_dVi;


726:   MatZeroEntries(B);
727:   DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
728:   DMCompositeScatter(user->dmpgrid,X,Xgen,Xnet);

730:   VecGetArray(Xgen,&xgen);
731:   VecGetArray(Xnet,&xnet);

733:   /* Generator subsystem */
734:   for (i=0; i < ngen; i++) {
735:     Eqp   = xgen[idx];
736:     Edp   = xgen[idx+1];
737:     delta = xgen[idx+2];
738:     Id    = xgen[idx+4];
739:     Iq    = xgen[idx+5];
740:     Efd   = xgen[idx+6];

742:     /*    fgen[idx]   = (-Eqp - (Xd[i] - Xdp[i])*Id + Efd)/Td0p[i]; */
743:     row[0] = idx;
744:     col[0] = idx;           col[1] = idx+4;          col[2] = idx+6;
745:     val[0] = -1/ Td0p[i]; val[1] = -(Xd[i] - Xdp[i])/ Td0p[i]; val[2] = 1/Td0p[i];

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

749:     /*    fgen[idx+1] = (-Edp + (Xq[i] - Xqp[i])*Iq)/Tq0p[i]; */
750:     row[0] = idx + 1;
751:     col[0] = idx + 1;       col[1] = idx+5;
752:     val[0] = -1/Tq0p[i]; val[1] = (Xq[i] - Xqp[i])/Tq0p[i];
753:     MatSetValues(J,1,row,2,col,val,INSERT_VALUES);

755:     /*    fgen[idx+2] = w - w_s; */
756:     row[0] = idx + 2;
757:     col[0] = idx + 2; col[1] = idx + 3;
758:     val[0] = 0;       val[1] = 1;
759:     MatSetValues(J,1,row,2,col,val,INSERT_VALUES);

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

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

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

773:     Zdq_inv[0] = Rs[i]/det;
774:     Zdq_inv[1] = Xqp[i]/det;
775:     Zdq_inv[2] = -Xdp[i]/det;
776:     Zdq_inv[3] = Rs[i]/det;

778:     dVd_dVr    = PetscSinScalar(delta); dVd_dVi = -PetscCosScalar(delta);
779:     dVq_dVr    = PetscCosScalar(delta); dVq_dVi = PetscSinScalar(delta);
780:     dVd_ddelta = Vr*PetscCosScalar(delta) + Vi*PetscSinScalar(delta);
781:     dVq_ddelta = -Vr*PetscSinScalar(delta) + Vi*PetscCosScalar(delta);

783:     /*    fgen[idx+4] = Zdq_inv[0]*(-Edp + Vd) + Zdq_inv[1]*(-Eqp + Vq) + Id; */
784:     row[0] = idx+4;
785:     col[0] = idx;         col[1] = idx+1;        col[2] = idx + 2;
786:     val[0] = -Zdq_inv[1]; val[1] = -Zdq_inv[0];  val[2] = Zdq_inv[0]*dVd_ddelta + Zdq_inv[1]*dVq_ddelta;
787:     col[3] = idx + 4; col[4] = net_start+2*gbus[i];                     col[5] = net_start + 2*gbus[i]+1;
788:     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;
789:     MatSetValues(J,1,row,6,col,val,INSERT_VALUES);

791:     /*  fgen[idx+5] = Zdq_inv[2]*(-Edp + Vd) + Zdq_inv[3]*(-Eqp + Vq) + Iq; */
792:     row[0] = idx+5;
793:     col[0] = idx;         col[1] = idx+1;        col[2] = idx + 2;
794:     val[0] = -Zdq_inv[3]; val[1] = -Zdq_inv[2];  val[2] = Zdq_inv[2]*dVd_ddelta + Zdq_inv[3]*dVq_ddelta;
795:     col[3] = idx + 5; col[4] = net_start+2*gbus[i];                     col[5] = net_start + 2*gbus[i]+1;
796:     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;
797:     MatSetValues(J,1,row,6,col,val,INSERT_VALUES);

799:     dIGr_ddelta = Id*PetscCosScalar(delta) - Iq*PetscSinScalar(delta);
800:     dIGi_ddelta = Id*PetscSinScalar(delta) + Iq*PetscCosScalar(delta);
801:     dIGr_dId    = PetscSinScalar(delta);  dIGr_dIq = PetscCosScalar(delta);
802:     dIGi_dId    = -PetscCosScalar(delta); dIGi_dIq = PetscSinScalar(delta);

804:     /* fnet[2*gbus[i]]   -= IGi; */
805:     row[0] = net_start + 2*gbus[i];
806:     col[0] = idx+2;        col[1] = idx + 4;   col[2] = idx + 5;
807:     val[0] = -dIGi_ddelta; val[1] = -dIGi_dId; val[2] = -dIGi_dIq;
808:     MatSetValues(J,1,row,3,col,val,INSERT_VALUES);

810:     /* fnet[2*gbus[i]+1]   -= IGr; */
811:     row[0] = net_start + 2*gbus[i]+1;
812:     col[0] = idx+2;        col[1] = idx + 4;   col[2] = idx + 5;
813:     val[0] = -dIGr_ddelta; val[1] = -dIGr_dId; val[2] = -dIGr_dIq;
814:     MatSetValues(J,1,row,3,col,val,INSERT_VALUES);

816:     Vm = PetscSqrtScalar(Vd*Vd + Vq*Vq);

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

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

823:     row[0] = idx + 6;
824:     col[0] = idx + 6;                     col[1] = idx + 8;
825:     val[0] = (-KE[i] - dSE_dEfd)/TE[i];  val[1] = 1/TE[i];
826:     MatSetValues(J,1,row,2,col,val,INSERT_VALUES);

828:     /* Exciter differential equations */

830:     /*    fgen[idx+7] = (-RF + KF[i]*Efd/TF[i])/TF[i]; */
831:     row[0] = idx + 7;
832:     col[0] = idx + 6;       col[1] = idx + 7;
833:     val[0] = (KF[i]/TF[i])/TF[i];  val[1] = -1/TF[i];
834:     MatSetValues(J,1,row,2,col,val,INSERT_VALUES);

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

839:     row[0] = idx + 8;
840:     if(VRatmax[i]) {
841:       col[0] = idx + 8; val[0] = 1.0;
842:       MatSetValues(J,1,row,1,col,val,INSERT_VALUES);
843:     } else if(VRatmin[i]) {
844:       col[0] = idx + 8; val[0] = -1.0;
845:       MatSetValues(J,1,row,1,col,val,INSERT_VALUES);
846:     } else {
847:       dVm_dVd    = Vd/Vm; dVm_dVq = Vq/Vm;
848:       dVm_dVr    = dVm_dVd*dVd_dVr + dVm_dVq*dVq_dVr;
849:       dVm_dVi    = dVm_dVd*dVd_dVi + dVm_dVq*dVq_dVi;
850:       row[0]     = idx + 8;
851:       col[0]     = idx + 6;           col[1] = idx + 7; col[2] = idx + 8;
852:       val[0]     = -(KA[i]*KF[i]/TF[i])/TA[i]; val[1] = KA[i]/TA[i];  val[2] = -1/TA[i];
853:       col[3]     = net_start + 2*gbus[i]; col[4] = net_start + 2*gbus[i]+1;
854:       val[3]     = -KA[i]*dVm_dVr/TA[i];         val[4] = -KA[i]*dVm_dVi/TA[i];
855:       MatSetValues(J,1,row,5,col,val,INSERT_VALUES);
856:     }

858:     idx        = idx + 9;
859:   }

861:   for (i=0; i<nbus; i++) {
862:     MatGetRow(user->Ybus,2*i,&ncols,&cols,&yvals);
863:     row[0] = net_start + 2*i;
864:     for (k=0; k<ncols; k++) {
865:       col[k] = net_start + cols[k];
866:       val[k] = yvals[k];
867:     }
868:     MatSetValues(J,1,row,ncols,col,val,INSERT_VALUES);
869:     MatRestoreRow(user->Ybus,2*i,&ncols,&cols,&yvals);

871:     MatGetRow(user->Ybus,2*i+1,&ncols,&cols,&yvals);
872:     row[0] = net_start + 2*i+1;
873:     for (k=0; k<ncols; k++) {
874:       col[k] = net_start + cols[k];
875:       val[k] = yvals[k];
876:     }
877:     MatSetValues(J,1,row,ncols,col,val,INSERT_VALUES);
878:     MatRestoreRow(user->Ybus,2*i+1,&ncols,&cols,&yvals);
879:   }

881:   MatAssemblyBegin(J,MAT_FLUSH_ASSEMBLY);
882:   MatAssemblyEnd(J,MAT_FLUSH_ASSEMBLY);

884:   VecGetArray(user->V0,&v0);
885:   for (i=0; i < nload; i++) {
886:     Vr      = xnet[2*lbus[i]]; /* Real part of load bus voltage */
887:     Vi      = xnet[2*lbus[i]+1]; /* Imaginary part of the load bus voltage */
888:     Vm      = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm; Vm4 = Vm2*Vm2;
889:     Vm0     = PetscSqrtScalar(v0[2*lbus[i]]*v0[2*lbus[i]] + v0[2*lbus[i]+1]*v0[2*lbus[i]+1]);
890:     PD      = QD = 0.0;
891:     dPD_dVr = dPD_dVi = dQD_dVr = dQD_dVi = 0.0;
892:     for (k=0; k < ld_nsegsp[i]; k++) {
893:       PD      += ld_alphap[k]*PD0[i]*PetscPowScalar((Vm/Vm0),ld_betap[k]);
894:       dPD_dVr += ld_alphap[k]*ld_betap[k]*PD0[i]*PetscPowScalar((1/Vm0),ld_betap[k])*Vr*PetscPowScalar(Vm,(ld_betap[k]-2));
895:       dPD_dVi += ld_alphap[k]*ld_betap[k]*PD0[i]*PetscPowScalar((1/Vm0),ld_betap[k])*Vi*PetscPowScalar(Vm,(ld_betap[k]-2));
896:     }
897:     for (k=0; k < ld_nsegsq[i]; k++) {
898:       QD      += ld_alphaq[k]*QD0[i]*PetscPowScalar((Vm/Vm0),ld_betaq[k]);
899:       dQD_dVr += ld_alphaq[k]*ld_betaq[k]*QD0[i]*PetscPowScalar((1/Vm0),ld_betaq[k])*Vr*PetscPowScalar(Vm,(ld_betaq[k]-2));
900:       dQD_dVi += ld_alphaq[k]*ld_betaq[k]*QD0[i]*PetscPowScalar((1/Vm0),ld_betaq[k])*Vi*PetscPowScalar(Vm,(ld_betaq[k]-2));
901:     }

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

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

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


913:     /*    fnet[2*lbus[i]]   += IDi; */
914:     row[0] = net_start + 2*lbus[i];
915:     col[0] = net_start + 2*lbus[i];  col[1] = net_start + 2*lbus[i]+1;
916:     val[0] = dIDi_dVr;               val[1] = dIDi_dVi;
917:     MatSetValues(J,1,row,2,col,val,ADD_VALUES);
918:     /*    fnet[2*lbus[i]+1] += IDr; */
919:     row[0] = net_start + 2*lbus[i]+1;
920:     col[0] = net_start + 2*lbus[i];  col[1] = net_start + 2*lbus[i]+1;
921:     val[0] = dIDr_dVr;               val[1] = dIDr_dVi;
922:     MatSetValues(J,1,row,2,col,val,ADD_VALUES);
923:   }
924:   VecRestoreArray(user->V0,&v0);

926:   VecRestoreArray(Xgen,&xgen);
927:   VecRestoreArray(Xnet,&xnet);

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

931:   MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
932:   MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
933:   return(0);
934: }

936: /*
937:    J = [I, 0
938:         dg_dx, dg_dy]
939: */
940: PetscErrorCode AlgJacobian(SNES snes,Vec X,Mat A,Mat B,void *ctx)
941: {
943:   Userctx        *user=(Userctx*)ctx;

946:   ResidualJacobian(X,A,B,ctx);
947:   MatSetOption(A,MAT_KEEP_NONZERO_PATTERN,PETSC_TRUE);
948:   MatZeroRowsIS(A,user->is_diff,1.0,NULL,NULL);
949:   return(0);
950: }

952: /*
953:    J = [-df_dx, -df_dy
954:         dg_dx, dg_dy]
955: */

957: PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec X,Mat A,Mat B,void *ctx)
958: {
960:   Userctx        *user=(Userctx*)ctx;

963:   user->t = t;

965:   ResidualJacobian(X,A,B,user);

967:   return(0);
968: }

970: /*
971:    J = [df_dx-aI, df_dy
972:         dg_dx, dg_dy]
973: */

975: PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,Userctx *user)
976: {
978:   PetscScalar    atmp = (PetscScalar) a;
979:   PetscInt       i,row;

982:   user->t = t;
983:   atmp *= -1;

985:   RHSJacobian(ts,t,X,A,B,user);
986:   for (i=0;i < ngen;i++) {
987:     row = 9*i;
988:     MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
989:     row  = 9*i+1;
990:     MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
991:     row  = 9*i+2;
992:     MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
993:     row  = 9*i+3;
994:     MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
995:     row  = 9*i+6;
996:     MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
997:     row  = 9*i+7;
998:     MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
999:     row  = 9*i+8;
1000:     MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
1001:   }
1002:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
1003:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
1004:   return(0);
1005: }

1007: int main(int argc,char **argv)
1008: {
1009:   TS             ts;
1010:   SNES           snes_alg;
1012:   PetscMPIInt    size;
1013:   Userctx        user;
1014:   PetscViewer    Xview,Ybusview,viewer;
1015:   Vec            X,F_alg;
1016:   Mat            J,A;
1017:   PetscInt       i,idx,*idx2;
1018:   Vec            Xdot;
1019:   PetscScalar    *x,*mat,*amat;
1020:   Vec            vatol;
1021:   PetscInt       *direction;
1022:   PetscBool      *terminate;
1023:   const PetscInt *idx3;
1024:   PetscScalar    *vatoli;
1025:   PetscInt       k;


1028:   PetscInitialize(&argc,&argv,"petscoptions",help);
1029:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
1030:   if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs");

1032:   user.neqs_gen   = 9*ngen; /* # eqs. for generator subsystem */
1033:   user.neqs_net   = 2*nbus; /* # eqs. for network subsystem   */
1034:   user.neqs_pgrid = user.neqs_gen + user.neqs_net;

1036:   /* Create indices for differential and algebraic equations */

1038:   PetscMalloc1(7*ngen,&idx2);
1039:   for (i=0; i<ngen; i++) {
1040:     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;
1041:     idx2[7*i+4] = 9*i+6; idx2[7*i+5] = 9*i+7; idx2[7*i+6] = 9*i+8;
1042:   }
1043:   ISCreateGeneral(PETSC_COMM_WORLD,7*ngen,idx2,PETSC_COPY_VALUES,&user.is_diff);
1044:   ISComplement(user.is_diff,0,user.neqs_pgrid,&user.is_alg);
1045:   PetscFree(idx2);

1047:   /* Read initial voltage vector and Ybus */
1048:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"X.bin",FILE_MODE_READ,&Xview);
1049:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"Ybus.bin",FILE_MODE_READ,&Ybusview);

1051:   VecCreate(PETSC_COMM_WORLD,&user.V0);
1052:   VecSetSizes(user.V0,PETSC_DECIDE,user.neqs_net);
1053:   VecLoad(user.V0,Xview);

1055:   MatCreate(PETSC_COMM_WORLD,&user.Ybus);
1056:   MatSetSizes(user.Ybus,PETSC_DECIDE,PETSC_DECIDE,user.neqs_net,user.neqs_net);
1057:   MatSetType(user.Ybus,MATBAIJ);
1058:   /*  MatSetBlockSize(user.Ybus,2); */
1059:   MatLoad(user.Ybus,Ybusview);

1061:   /* Set run time options */
1062:   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Transient stability fault options","");
1063:   {
1064:     user.tfaulton  = 1.0;
1065:     user.tfaultoff = 1.2;
1066:     user.Rfault    = 0.0001;
1067:     user.setisdiff = PETSC_FALSE;
1068:     user.semiexplicit = PETSC_FALSE;
1069:     user.faultbus  = 8;
1070:     PetscOptionsReal("-tfaulton","","",user.tfaulton,&user.tfaulton,NULL);
1071:     PetscOptionsReal("-tfaultoff","","",user.tfaultoff,&user.tfaultoff,NULL);
1072:     PetscOptionsInt("-faultbus","","",user.faultbus,&user.faultbus,NULL);
1073:     user.t0        = 0.0;
1074:     user.tmax      = 5.0;
1075:     PetscOptionsReal("-t0","","",user.t0,&user.t0,NULL);
1076:     PetscOptionsReal("-tmax","","",user.tmax,&user.tmax,NULL);
1077:     PetscOptionsBool("-setisdiff","","",user.setisdiff,&user.setisdiff,NULL);
1078:     PetscOptionsBool("-dae_semiexplicit","","",user.semiexplicit,&user.semiexplicit,NULL);
1079:   }
1080:   PetscOptionsEnd();

1082:   PetscViewerDestroy(&Xview);
1083:   PetscViewerDestroy(&Ybusview);

1085:   /* Create DMs for generator and network subsystems */
1086:   DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,user.neqs_gen,1,1,NULL,&user.dmgen);
1087:   DMSetOptionsPrefix(user.dmgen,"dmgen_");
1088:   DMSetFromOptions(user.dmgen);
1089:   DMSetUp(user.dmgen);
1090:   DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,user.neqs_net,1,1,NULL,&user.dmnet);
1091:   DMSetOptionsPrefix(user.dmnet,"dmnet_");
1092:   DMSetFromOptions(user.dmnet);
1093:   DMSetUp(user.dmnet);
1094:   /* Create a composite DM packer and add the two DMs */
1095:   DMCompositeCreate(PETSC_COMM_WORLD,&user.dmpgrid);
1096:   DMSetOptionsPrefix(user.dmpgrid,"pgrid_");
1097:   DMCompositeAddDM(user.dmpgrid,user.dmgen);
1098:   DMCompositeAddDM(user.dmpgrid,user.dmnet);

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

1102:   MatCreate(PETSC_COMM_WORLD,&J);
1103:   MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,user.neqs_pgrid,user.neqs_pgrid);
1104:   MatSetFromOptions(J);
1105:   PreallocateJacobian(J,&user);

1107:   /* Create matrix to save solutions at each time step */
1108:   user.stepnum = 0;

1110:   MatCreateSeqDense(PETSC_COMM_SELF,user.neqs_pgrid+1,1002,NULL,&user.Sol);
1111:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1112:      Create timestepping solver context
1113:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1114:   TSCreate(PETSC_COMM_WORLD,&ts);
1115:   TSSetProblemType(ts,TS_NONLINEAR);
1116:   if(user.semiexplicit) {
1117:     TSSetType(ts,TSRK);
1118:     TSSetRHSFunction(ts,NULL,RHSFunction,&user);
1119:     TSSetRHSJacobian(ts,J,J,RHSJacobian,&user);
1120:   } else {
1121:     TSSetType(ts,TSCN);
1122:     TSSetEquationType(ts,TS_EQ_DAE_IMPLICIT_INDEX1);
1123:     TSARKIMEXSetFullyImplicit(ts,PETSC_TRUE);
1124:     TSSetIFunction(ts,NULL,(TSIFunction) IFunction,&user);
1125:     TSSetIJacobian(ts,J,J,(TSIJacobian)IJacobian,&user);
1126:   }
1127:   TSSetApplicationContext(ts,&user);

1129:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1130:      Set initial conditions
1131:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1132:   SetInitialGuess(X,&user);
1133:   /* Just to set up the Jacobian structure */

1135:   VecDuplicate(X,&Xdot);
1136:   IJacobian(ts,0.0,X,Xdot,0.0,J,J,&user);
1137:   VecDestroy(&Xdot);

1139:   /* Save initial solution */

1141:   idx=user.stepnum*(user.neqs_pgrid+1);
1142:   MatDenseGetArray(user.Sol,&mat);
1143:   VecGetArray(X,&x);

1145:   mat[idx] = 0.0;

1147:   PetscMemcpy(mat+idx+1,x,user.neqs_pgrid*sizeof(PetscScalar));
1148:   MatDenseRestoreArray(user.Sol,&mat);
1149:   VecRestoreArray(X,&x);
1150:   user.stepnum++;

1152:   TSSetMaxTime(ts,user.tmax);
1153:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);
1154:   TSSetTimeStep(ts,0.01);
1155:   TSSetFromOptions(ts);
1156:   TSSetPostStep(ts,SaveSolution);
1157:   TSSetSolution(ts,X);

1159:   PetscMalloc1((2*ngen+2),&direction);
1160:   PetscMalloc1((2*ngen+2),&terminate);
1161:   direction[0] = direction[1] = 1;
1162:   terminate[0] = terminate[1] = PETSC_FALSE;
1163:   for (i=0; i < ngen;i++) {
1164:     direction[2+2*i] = -1; direction[2+2*i+1] = 1;
1165:     terminate[2+2*i] = terminate[2+2*i+1] = PETSC_FALSE;
1166:   }

1168:   TSSetEventHandler(ts,2*ngen+2,direction,terminate,EventFunction,PostEventFunction,(void*)&user);

1170:   if(user.semiexplicit) {
1171:     /* Use a semi-explicit approach with the time-stepping done by an explicit method and the
1172:        algrebraic part solved via PostStage and PostEvaluate callbacks 
1173:     */
1174:     TSSetType(ts,TSRK);
1175:     TSSetPostStage(ts,PostStage);
1176:     TSSetPostEvaluate(ts,PostEvaluate);
1177:   }


1180:   if(user.setisdiff) {
1181:     /* Create vector of absolute tolerances and set the algebraic part to infinity */
1182:     VecDuplicate(X,&vatol);
1183:     VecSet(vatol,100000.0);
1184:     VecGetArray(vatol,&vatoli);
1185:     ISGetIndices(user.is_diff,&idx3);
1186:     for(k=0; k < 7*ngen; k++) vatoli[idx3[k]] = 1e-2;
1187:     VecRestoreArray(vatol,&vatoli);
1188:   }

1190:   /* Create the nonlinear solver for solving the algebraic system */
1191:   /* Note that although the algebraic system needs to be solved only for
1192:      Idq and V, we reuse the entire system including xgen. The xgen
1193:      variables are held constant by setting their residuals to 0 and
1194:      putting a 1 on the Jacobian diagonal for xgen rows
1195:   */

1197:   VecDuplicate(X,&F_alg);
1198:   SNESCreate(PETSC_COMM_WORLD,&snes_alg);
1199:   SNESSetFunction(snes_alg,F_alg,AlgFunction,&user);
1200:   SNESSetJacobian(snes_alg,J,J,AlgJacobian,&user);

1202:   SNESSetFromOptions(snes_alg);

1204:   user.snes_alg=snes_alg;
1205: 
1206:   /* Solve */
1207:   TSSolve(ts,X);

1209:   MatAssemblyBegin(user.Sol,MAT_FINAL_ASSEMBLY);
1210:   MatAssemblyEnd(user.Sol,MAT_FINAL_ASSEMBLY);

1212:   MatCreateSeqDense(PETSC_COMM_SELF,user.neqs_pgrid+1,user.stepnum,NULL,&A);
1213:   MatDenseGetArray(user.Sol,&mat);
1214:   MatDenseGetArray(A,&amat);
1215:   PetscMemcpy(amat,mat,(user.stepnum*(user.neqs_pgrid+1))*sizeof(PetscScalar));
1216:   MatDenseRestoreArray(A,&amat);
1217:   MatDenseRestoreArray(user.Sol,&mat);
1218:   PetscViewerBinaryOpen(PETSC_COMM_SELF,"out.bin",FILE_MODE_WRITE,&viewer);
1219:   MatView(A,viewer);
1220:   PetscViewerDestroy(&viewer);
1221:   MatDestroy(&A);

1223:   PetscFree(direction);
1224:   PetscFree(terminate);
1225:   SNESDestroy(&snes_alg);
1226:   VecDestroy(&F_alg);
1227:   MatDestroy(&J);
1228:   MatDestroy(&user.Ybus);
1229:   MatDestroy(&user.Sol);
1230:   VecDestroy(&X);
1231:   VecDestroy(&user.V0);
1232:   DMDestroy(&user.dmgen);
1233:   DMDestroy(&user.dmnet);
1234:   DMDestroy(&user.dmpgrid);
1235:   ISDestroy(&user.is_diff);
1236:   ISDestroy(&user.is_alg);
1237:   TSDestroy(&ts);
1238:   if(user.setisdiff) {
1239:     VecDestroy(&vatol);
1240:   }
1241:   PetscFinalize();
1242:   return ierr;
1243: }