Actual source code: ex9bus.c


  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 coordinates.\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: */

 42: #include <petscts.h>
 43: #include <petscdm.h>
 44: #include <petscdmda.h>
 45: #include <petscdmcomposite.h>

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

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

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

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

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

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

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

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


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

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

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

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

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

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

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

183:   return 0;
184: }

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


199:   DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
200:   DMCompositeScatter(user->dmpgrid,X,Xgen,Xnet);

202:   VecGetArray(Xgen,&xgen);
203:   VecGetArray(Xnet,&xnet);

205:   for (i=0; i < nevents; i++) {
206:     if (event_list[i] == 0) {
207:       /* Apply disturbance - resistive fault at user->faultbus */
208:       /* This is done by adding shunt conductance to the diagonal location
209:          in the Ybus matrix */
210:       row_loc = 2*user->faultbus; col_loc = 2*user->faultbus+1; /* Location for G */
211:       val     = 1/user->Rfault;
212:       MatSetValues(user->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);
213:       row_loc = 2*user->faultbus+1; col_loc = 2*user->faultbus; /* Location for G */
214:       val     = 1/user->Rfault;
215:       MatSetValues(user->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);

217:       MatAssemblyBegin(user->Ybus,MAT_FINAL_ASSEMBLY);
218:       MatAssemblyEnd(user->Ybus,MAT_FINAL_ASSEMBLY);

220:       /* Solve the algebraic equations */
221:       SNESSolve(user->snes_alg,NULL,X);
222:     } else if (event_list[i] == 1) {
223:       /* Remove the fault */
224:       row_loc = 2*user->faultbus; col_loc = 2*user->faultbus+1;
225:       val     = -1/user->Rfault;
226:       MatSetValues(user->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);
227:       row_loc = 2*user->faultbus+1; col_loc = 2*user->faultbus;
228:       val     = -1/user->Rfault;
229:       MatSetValues(user->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);

231:       MatAssemblyBegin(user->Ybus,MAT_FINAL_ASSEMBLY);
232:       MatAssemblyEnd(user->Ybus,MAT_FINAL_ASSEMBLY);

234:       /* Solve the algebraic equations */
235:       SNESSolve(user->snes_alg,NULL,X);

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

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

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

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

289:   VecRestoreArray(Xgen,&xgen);
290:   VecRestoreArray(Xnet,&xnet);

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

294:   return 0;
295: }

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

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

313: /* Saves the solution at each time to a matrix */
314: PetscErrorCode SaveSolution(TS ts)
315: {
316:   Userctx           *user;
317:   Vec               X;
318:   const PetscScalar *x;
319:   PetscScalar       *mat;
320:   PetscInt          idx;
321:   PetscReal         t;

323:   TSGetApplicationContext(ts,&user);
324:   TSGetTime(ts,&t);
325:   TSGetSolution(ts,&X);
326:   idx      = user->stepnum*(user->neqs_pgrid+1);
327:   MatDenseGetArray(user->Sol,&mat);
328:   VecGetArrayRead(X,&x);
329:   mat[idx] = t;
330:   PetscArraycpy(mat+idx+1,x,user->neqs_pgrid);
331:   MatDenseRestoreArray(user->Sol,&mat);
332:   VecRestoreArrayRead(X,&x);
333:   user->stepnum++;
334:   return 0;
335: }

337: PetscErrorCode SetInitialGuess(Vec X,Userctx *user)
338: {
339:   Vec               Xgen,Xnet;
340:   PetscScalar       *xgen;
341:   const PetscScalar *xnet;
342:   PetscInt          i,idx=0;
343:   PetscScalar       Vr,Vi,IGr,IGi,Vm,Vm2;
344:   PetscScalar       Eqp,Edp,delta;
345:   PetscScalar       Efd,RF,VR; /* Exciter variables */
346:   PetscScalar       Id,Iq;  /* Generator dq axis currents */
347:   PetscScalar       theta,Vd,Vq,SE;

349:   M[0] = 2*H[0]/w_s; M[1] = 2*H[1]/w_s; M[2] = 2*H[2]/w_s;
350:   D[0] = 0.1*M[0]; D[1] = 0.1*M[1]; D[2] = 0.1*M[2];

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

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

357:   /* Generator subsystem initialization */
358:   VecGetArrayWrite(Xgen,&xgen);
359:   VecGetArrayRead(Xnet,&xnet);

361:   for (i=0; i < ngen; i++) {
362:     Vr  = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */
363:     Vi  = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */
364:     Vm  = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm;
365:     IGr = (Vr*PG[i] + Vi*QG[i])/Vm2;
366:     IGi = (Vi*PG[i] - Vr*QG[i])/Vm2;

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

370:     theta = PETSC_PI/2.0 - delta;

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

375:     Vd = Vr*PetscCosScalar(theta) - Vi*PetscSinScalar(theta);
376:     Vq = Vr*PetscSinScalar(theta) + Vi*PetscCosScalar(theta);

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

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

383:     /* The generator variables are ordered as [Eqp,Edp,delta,w,Id,Iq] */
384:     xgen[idx]   = Eqp;
385:     xgen[idx+1] = Edp;
386:     xgen[idx+2] = delta;
387:     xgen[idx+3] = w_s;

389:     idx = idx + 4;

391:     xgen[idx]   = Id;
392:     xgen[idx+1] = Iq;

394:     idx = idx + 2;

396:     /* Exciter */
397:     Efd = Eqp + (Xd[i] - Xdp[i])*Id;
398:     SE  = k1[i]*PetscExpScalar(k2[i]*Efd);
399:     VR  =  KE[i]*Efd + SE;
400:     RF  =  KF[i]*Efd/TF[i];

402:     xgen[idx]   = Efd;
403:     xgen[idx+1] = RF;
404:     xgen[idx+2] = VR;

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

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

410:     idx = idx + 3;
411:   }
412:   VecRestoreArrayWrite(Xgen,&xgen);
413:   VecRestoreArrayRead(Xnet,&xnet);

415:   /* VecView(Xgen,0); */
416:   DMCompositeGather(user->dmpgrid,INSERT_VALUES,X,Xgen,Xnet);
417:   DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);
418:   return 0;
419: }

421: /* Computes F = [f(x,y);g(x,y)] */
422: PetscErrorCode ResidualFunction(Vec X, Vec F, Userctx *user)
423: {
424:   Vec               Xgen,Xnet,Fgen,Fnet;
425:   const PetscScalar *xgen,*xnet;
426:   PetscScalar       *fgen,*fnet;
427:   PetscInt          i,idx=0;
428:   PetscScalar       Vr,Vi,Vm,Vm2;
429:   PetscScalar       Eqp,Edp,delta,w; /* Generator variables */
430:   PetscScalar       Efd,RF,VR; /* Exciter variables */
431:   PetscScalar       Id,Iq;  /* Generator dq axis currents */
432:   PetscScalar       Vd,Vq,SE;
433:   PetscScalar       IGr,IGi,IDr,IDi;
434:   PetscScalar       Zdq_inv[4],det;
435:   PetscScalar       PD,QD,Vm0,*v0;
436:   PetscInt          k;

438:   VecZeroEntries(F);
439:   DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
440:   DMCompositeGetLocalVectors(user->dmpgrid,&Fgen,&Fnet);
441:   DMCompositeScatter(user->dmpgrid,X,Xgen,Xnet);
442:   DMCompositeScatter(user->dmpgrid,F,Fgen,Fnet);

444:   /* Network current balance residual IG + Y*V + IL = 0. Only YV is added here.
445:      The generator current injection, IG, and load current injection, ID are added later
446:   */
447:   /* Note that the values in Ybus are stored assuming the imaginary current balance
448:      equation is ordered first followed by real current balance equation for each bus.
449:      Thus imaginary current contribution goes in location 2*i, and
450:      real current contribution in 2*i+1
451:   */
452:   MatMult(user->Ybus,Xnet,Fnet);

454:   VecGetArrayRead(Xgen,&xgen);
455:   VecGetArrayRead(Xnet,&xnet);
456:   VecGetArrayWrite(Fgen,&fgen);
457:   VecGetArrayWrite(Fnet,&fnet);

459:   /* Generator subsystem */
460:   for (i=0; i < ngen; i++) {
461:     Eqp   = xgen[idx];
462:     Edp   = xgen[idx+1];
463:     delta = xgen[idx+2];
464:     w     = xgen[idx+3];
465:     Id    = xgen[idx+4];
466:     Iq    = xgen[idx+5];
467:     Efd   = xgen[idx+6];
468:     RF    = xgen[idx+7];
469:     VR    = xgen[idx+8];

471:     /* Generator differential equations */
472:     fgen[idx]   = (-Eqp - (Xd[i] - Xdp[i])*Id + Efd)/Td0p[i];
473:     fgen[idx+1] = (-Edp + (Xq[i] - Xqp[i])*Iq)/Tq0p[i];
474:     fgen[idx+2] = w - w_s;
475:     fgen[idx+3] = (TM[i] - Edp*Id - Eqp*Iq - (Xqp[i] - Xdp[i])*Id*Iq - D[i]*(w - w_s))/M[i];

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

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

484:     Zdq_inv[0] = Rs[i]/det;
485:     Zdq_inv[1] = Xqp[i]/det;
486:     Zdq_inv[2] = -Xdp[i]/det;
487:     Zdq_inv[3] = Rs[i]/det;

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

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

495:     fnet[2*gbus[i]]   -= IGi;
496:     fnet[2*gbus[i]+1] -= IGr;

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

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

502:     /* Exciter differential equations */
503:     fgen[idx+6] = (-KE[i]*Efd - SE + VR)/TE[i];
504:     fgen[idx+7] = (-RF + KF[i]*Efd/TF[i])/TF[i];
505:     if (VRatmax[i]) fgen[idx+8] = VR - VRMAX[i];
506:     else if (VRatmin[i]) fgen[idx+8] = VRMIN[i] - VR;
507:     else fgen[idx+8] = (-VR + KA[i]*RF - KA[i]*KF[i]*Efd/TF[i] + KA[i]*(Vref[i] - Vm))/TA[i];

509:     idx = idx + 9;
510:   }

512:   VecGetArray(user->V0,&v0);
513:   for (i=0; i < nload; i++) {
514:     Vr  = xnet[2*lbus[i]]; /* Real part of load bus voltage */
515:     Vi  = xnet[2*lbus[i]+1]; /* Imaginary part of the load bus voltage */
516:     Vm  = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm;
517:     Vm0 = PetscSqrtScalar(v0[2*lbus[i]]*v0[2*lbus[i]] + v0[2*lbus[i]+1]*v0[2*lbus[i]+1]);
518:     PD  = QD = 0.0;
519:     for (k=0; k < ld_nsegsp[i]; k++) PD += ld_alphap[k]*PD0[i]*PetscPowScalar((Vm/Vm0),ld_betap[k]);
520:     for (k=0; k < ld_nsegsq[i]; k++) QD += ld_alphaq[k]*QD0[i]*PetscPowScalar((Vm/Vm0),ld_betaq[k]);

522:     /* Load currents */
523:     IDr = (PD*Vr + QD*Vi)/Vm2;
524:     IDi = (-QD*Vr + PD*Vi)/Vm2;

526:     fnet[2*lbus[i]]   += IDi;
527:     fnet[2*lbus[i]+1] += IDr;
528:   }
529:   VecRestoreArray(user->V0,&v0);

531:   VecRestoreArrayRead(Xgen,&xgen);
532:   VecRestoreArrayRead(Xnet,&xnet);
533:   VecRestoreArrayWrite(Fgen,&fgen);
534:   VecRestoreArrayWrite(Fnet,&fnet);

536:   DMCompositeGather(user->dmpgrid,INSERT_VALUES,F,Fgen,Fnet);
537:   DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);
538:   DMCompositeRestoreLocalVectors(user->dmpgrid,&Fgen,&Fnet);
539:   return 0;
540: }

542: /*   f(x,y)
543:      g(x,y)
544:  */
545: PetscErrorCode RHSFunction(TS ts,PetscReal t, Vec X, Vec F, void *ctx)
546: {
547:   Userctx        *user=(Userctx*)ctx;

549:   user->t = t;
550:   ResidualFunction(X,F,user);
551:   return 0;
552: }

554: /* f(x,y) - \dot{x}
555:      g(x,y) = 0
556:  */
557: PetscErrorCode IFunction(TS ts,PetscReal t, Vec X, Vec Xdot, Vec F, void *ctx)
558: {
559:   PetscScalar       *f;
560:   const PetscScalar *xdot;
561:   PetscInt          i;

563:   RHSFunction(ts,t,X,F,ctx);
564:   VecScale(F,-1.0);
565:   VecGetArray(F,&f);
566:   VecGetArrayRead(Xdot,&xdot);
567:   for (i=0;i < ngen;i++) {
568:     f[9*i]   += xdot[9*i];
569:     f[9*i+1] += xdot[9*i+1];
570:     f[9*i+2] += xdot[9*i+2];
571:     f[9*i+3] += xdot[9*i+3];
572:     f[9*i+6] += xdot[9*i+6];
573:     f[9*i+7] += xdot[9*i+7];
574:     f[9*i+8] += xdot[9*i+8];
575:   }
576:   VecRestoreArray(F,&f);
577:   VecRestoreArrayRead(Xdot,&xdot);
578:   return 0;
579: }

581: /* This function is used for solving the algebraic system only during fault on and
582:    off times. It computes the entire F and then zeros out the part corresponding to
583:    differential equations
584:  F = [0;g(y)];
585: */
586: PetscErrorCode AlgFunction(SNES snes, Vec X, Vec F, void *ctx)
587: {
588:   Userctx        *user=(Userctx*)ctx;
589:   PetscScalar    *f;
590:   PetscInt       i;

592:   ResidualFunction(X,F,user);
593:   VecGetArray(F,&f);
594:   for (i=0; i < ngen; i++) {
595:     f[9*i]   = 0;
596:     f[9*i+1] = 0;
597:     f[9*i+2] = 0;
598:     f[9*i+3] = 0;
599:     f[9*i+6] = 0;
600:     f[9*i+7] = 0;
601:     f[9*i+8] = 0;
602:   }
603:   VecRestoreArray(F,&f);
604:   return 0;
605: }

607: PetscErrorCode PostStage(TS ts, PetscReal t, PetscInt i, Vec *X)
608: {
609:   Userctx        *user;

611:   TSGetApplicationContext(ts,&user);
612:   SNESSolve(user->snes_alg,NULL,X[i]);
613:   return 0;
614: }

616: PetscErrorCode PostEvaluate(TS ts)
617: {
618:   Userctx        *user;
619:   Vec            X;

621:   TSGetApplicationContext(ts,&user);
622:   TSGetSolution(ts,&X);
623:   SNESSolve(user->snes_alg,NULL,X);
624:   return 0;
625: }

627: PetscErrorCode PreallocateJacobian(Mat J, Userctx *user)
628: {
629:   PetscInt       *d_nnz;
630:   PetscInt       i,idx=0,start=0;
631:   PetscInt       ncols;

633:   PetscMalloc1(user->neqs_pgrid,&d_nnz);
634:   for (i=0; i<user->neqs_pgrid; i++) d_nnz[i] = 0;
635:   /* Generator subsystem */
636:   for (i=0; i < ngen; i++) {

638:     d_nnz[idx]   += 3;
639:     d_nnz[idx+1] += 2;
640:     d_nnz[idx+2] += 2;
641:     d_nnz[idx+3] += 5;
642:     d_nnz[idx+4] += 6;
643:     d_nnz[idx+5] += 6;

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

648:     d_nnz[idx+6] += 2;
649:     d_nnz[idx+7] += 2;
650:     d_nnz[idx+8] += 5;

652:     idx = idx + 9;
653:   }
654:   start = user->neqs_gen;

656:   for (i=0; i < nbus; i++) {
657:     MatGetRow(user->Ybus,2*i,&ncols,NULL,NULL);
658:     d_nnz[start+2*i]   += ncols;
659:     d_nnz[start+2*i+1] += ncols;
660:     MatRestoreRow(user->Ybus,2*i,&ncols,NULL,NULL);
661:   }
662:   MatSeqAIJSetPreallocation(J,0,d_nnz);
663:   PetscFree(d_nnz);
664:   return 0;
665: }

667: /*
668:    J = [df_dx, df_dy
669:         dg_dx, dg_dy]
670: */
671: PetscErrorCode ResidualJacobian(Vec X,Mat J,Mat B,void *ctx)
672: {
673:   Userctx           *user = (Userctx*)ctx;
674:   Vec               Xgen,Xnet;
675:   const PetscScalar *xgen,*xnet;
676:   PetscInt          i,idx=0;
677:   PetscScalar       Vr,Vi,Vm,Vm2;
678:   PetscScalar       Eqp,Edp,delta; /* Generator variables */
679:   PetscScalar       Efd;
680:   PetscScalar       Id,Iq;  /* Generator dq axis currents */
681:   PetscScalar       Vd,Vq;
682:   PetscScalar       val[10];
683:   PetscInt          row[2],col[10];
684:   PetscInt          net_start = user->neqs_gen;
685:   PetscScalar       Zdq_inv[4],det;
686:   PetscScalar       dVd_dVr,dVd_dVi,dVq_dVr,dVq_dVi,dVd_ddelta,dVq_ddelta;
687:   PetscScalar       dIGr_ddelta,dIGi_ddelta,dIGr_dId,dIGr_dIq,dIGi_dId,dIGi_dIq;
688:   PetscScalar       dSE_dEfd;
689:   PetscScalar       dVm_dVd,dVm_dVq,dVm_dVr,dVm_dVi;
690:   PetscInt          ncols;
691:   const PetscInt    *cols;
692:   const PetscScalar *yvals;
693:   PetscInt          k;
694:   PetscScalar       PD,QD,Vm0,Vm4;
695:   const PetscScalar *v0;
696:   PetscScalar       dPD_dVr,dPD_dVi,dQD_dVr,dQD_dVi;
697:   PetscScalar       dIDr_dVr,dIDr_dVi,dIDi_dVr,dIDi_dVi;

699:   MatZeroEntries(B);
700:   DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
701:   DMCompositeScatter(user->dmpgrid,X,Xgen,Xnet);

703:   VecGetArrayRead(Xgen,&xgen);
704:   VecGetArrayRead(Xnet,&xnet);

706:   /* Generator subsystem */
707:   for (i=0; i < ngen; i++) {
708:     Eqp   = xgen[idx];
709:     Edp   = xgen[idx+1];
710:     delta = xgen[idx+2];
711:     Id    = xgen[idx+4];
712:     Iq    = xgen[idx+5];
713:     Efd   = xgen[idx+6];

715:     /*    fgen[idx]   = (-Eqp - (Xd[i] - Xdp[i])*Id + Efd)/Td0p[i]; */
716:     row[0] = idx;
717:     col[0] = idx;           col[1] = idx+4;          col[2] = idx+6;
718:     val[0] = -1/ Td0p[i]; val[1] = -(Xd[i] - Xdp[i])/ Td0p[i]; val[2] = 1/Td0p[i];

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

722:     /*    fgen[idx+1] = (-Edp + (Xq[i] - Xqp[i])*Iq)/Tq0p[i]; */
723:     row[0] = idx + 1;
724:     col[0] = idx + 1;       col[1] = idx+5;
725:     val[0] = -1/Tq0p[i]; val[1] = (Xq[i] - Xqp[i])/Tq0p[i];
726:     MatSetValues(J,1,row,2,col,val,INSERT_VALUES);

728:     /*    fgen[idx+2] = w - w_s; */
729:     row[0] = idx + 2;
730:     col[0] = idx + 2; col[1] = idx + 3;
731:     val[0] = 0;       val[1] = 1;
732:     MatSetValues(J,1,row,2,col,val,INSERT_VALUES);

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

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

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

746:     Zdq_inv[0] = Rs[i]/det;
747:     Zdq_inv[1] = Xqp[i]/det;
748:     Zdq_inv[2] = -Xdp[i]/det;
749:     Zdq_inv[3] = Rs[i]/det;

751:     dVd_dVr    = PetscSinScalar(delta); dVd_dVi = -PetscCosScalar(delta);
752:     dVq_dVr    = PetscCosScalar(delta); dVq_dVi = PetscSinScalar(delta);
753:     dVd_ddelta = Vr*PetscCosScalar(delta) + Vi*PetscSinScalar(delta);
754:     dVq_ddelta = -Vr*PetscSinScalar(delta) + Vi*PetscCosScalar(delta);

756:     /*    fgen[idx+4] = Zdq_inv[0]*(-Edp + Vd) + Zdq_inv[1]*(-Eqp + Vq) + Id; */
757:     row[0] = idx+4;
758:     col[0] = idx;         col[1] = idx+1;        col[2] = idx + 2;
759:     val[0] = -Zdq_inv[1]; val[1] = -Zdq_inv[0];  val[2] = Zdq_inv[0]*dVd_ddelta + Zdq_inv[1]*dVq_ddelta;
760:     col[3] = idx + 4; col[4] = net_start+2*gbus[i];                     col[5] = net_start + 2*gbus[i]+1;
761:     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;
762:     MatSetValues(J,1,row,6,col,val,INSERT_VALUES);

764:     /*  fgen[idx+5] = Zdq_inv[2]*(-Edp + Vd) + Zdq_inv[3]*(-Eqp + Vq) + Iq; */
765:     row[0] = idx+5;
766:     col[0] = idx;         col[1] = idx+1;        col[2] = idx + 2;
767:     val[0] = -Zdq_inv[3]; val[1] = -Zdq_inv[2];  val[2] = Zdq_inv[2]*dVd_ddelta + Zdq_inv[3]*dVq_ddelta;
768:     col[3] = idx + 5; col[4] = net_start+2*gbus[i];                     col[5] = net_start + 2*gbus[i]+1;
769:     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;
770:     MatSetValues(J,1,row,6,col,val,INSERT_VALUES);

772:     dIGr_ddelta = Id*PetscCosScalar(delta) - Iq*PetscSinScalar(delta);
773:     dIGi_ddelta = Id*PetscSinScalar(delta) + Iq*PetscCosScalar(delta);
774:     dIGr_dId    = PetscSinScalar(delta);  dIGr_dIq = PetscCosScalar(delta);
775:     dIGi_dId    = -PetscCosScalar(delta); dIGi_dIq = PetscSinScalar(delta);

777:     /* fnet[2*gbus[i]]   -= IGi; */
778:     row[0] = net_start + 2*gbus[i];
779:     col[0] = idx+2;        col[1] = idx + 4;   col[2] = idx + 5;
780:     val[0] = -dIGi_ddelta; val[1] = -dIGi_dId; val[2] = -dIGi_dIq;
781:     MatSetValues(J,1,row,3,col,val,INSERT_VALUES);

783:     /* fnet[2*gbus[i]+1]   -= IGr; */
784:     row[0] = net_start + 2*gbus[i]+1;
785:     col[0] = idx+2;        col[1] = idx + 4;   col[2] = idx + 5;
786:     val[0] = -dIGr_ddelta; val[1] = -dIGr_dId; val[2] = -dIGr_dIq;
787:     MatSetValues(J,1,row,3,col,val,INSERT_VALUES);

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

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

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

796:     row[0] = idx + 6;
797:     col[0] = idx + 6;                     col[1] = idx + 8;
798:     val[0] = (-KE[i] - dSE_dEfd)/TE[i];  val[1] = 1/TE[i];
799:     MatSetValues(J,1,row,2,col,val,INSERT_VALUES);

801:     /* Exciter differential equations */

803:     /*    fgen[idx+7] = (-RF + KF[i]*Efd/TF[i])/TF[i]; */
804:     row[0] = idx + 7;
805:     col[0] = idx + 6;       col[1] = idx + 7;
806:     val[0] = (KF[i]/TF[i])/TF[i];  val[1] = -1/TF[i];
807:     MatSetValues(J,1,row,2,col,val,INSERT_VALUES);

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

812:     row[0] = idx + 8;
813:     if (VRatmax[i]) {
814:       col[0] = idx + 8; val[0] = 1.0;
815:       MatSetValues(J,1,row,1,col,val,INSERT_VALUES);
816:     } else if (VRatmin[i]) {
817:       col[0] = idx + 8; val[0] = -1.0;
818:       MatSetValues(J,1,row,1,col,val,INSERT_VALUES);
819:     } else {
820:       dVm_dVd    = Vd/Vm; dVm_dVq = Vq/Vm;
821:       dVm_dVr    = dVm_dVd*dVd_dVr + dVm_dVq*dVq_dVr;
822:       dVm_dVi    = dVm_dVd*dVd_dVi + dVm_dVq*dVq_dVi;
823:       row[0]     = idx + 8;
824:       col[0]     = idx + 6;           col[1] = idx + 7; col[2] = idx + 8;
825:       val[0]     = -(KA[i]*KF[i]/TF[i])/TA[i]; val[1] = KA[i]/TA[i];  val[2] = -1/TA[i];
826:       col[3]     = net_start + 2*gbus[i]; col[4] = net_start + 2*gbus[i]+1;
827:       val[3]     = -KA[i]*dVm_dVr/TA[i];         val[4] = -KA[i]*dVm_dVi/TA[i];
828:       MatSetValues(J,1,row,5,col,val,INSERT_VALUES);
829:     }
830:     idx        = idx + 9;
831:   }

833:   for (i=0; i<nbus; i++) {
834:     MatGetRow(user->Ybus,2*i,&ncols,&cols,&yvals);
835:     row[0] = net_start + 2*i;
836:     for (k=0; k<ncols; k++) {
837:       col[k] = net_start + cols[k];
838:       val[k] = yvals[k];
839:     }
840:     MatSetValues(J,1,row,ncols,col,val,INSERT_VALUES);
841:     MatRestoreRow(user->Ybus,2*i,&ncols,&cols,&yvals);

843:     MatGetRow(user->Ybus,2*i+1,&ncols,&cols,&yvals);
844:     row[0] = net_start + 2*i+1;
845:     for (k=0; k<ncols; k++) {
846:       col[k] = net_start + cols[k];
847:       val[k] = yvals[k];
848:     }
849:     MatSetValues(J,1,row,ncols,col,val,INSERT_VALUES);
850:     MatRestoreRow(user->Ybus,2*i+1,&ncols,&cols,&yvals);
851:   }

853:   MatAssemblyBegin(J,MAT_FLUSH_ASSEMBLY);
854:   MatAssemblyEnd(J,MAT_FLUSH_ASSEMBLY);

856:   VecGetArrayRead(user->V0,&v0);
857:   for (i=0; i < nload; i++) {
858:     Vr      = xnet[2*lbus[i]]; /* Real part of load bus voltage */
859:     Vi      = xnet[2*lbus[i]+1]; /* Imaginary part of the load bus voltage */
860:     Vm      = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm; Vm4 = Vm2*Vm2;
861:     Vm0     = PetscSqrtScalar(v0[2*lbus[i]]*v0[2*lbus[i]] + v0[2*lbus[i]+1]*v0[2*lbus[i]+1]);
862:     PD      = QD = 0.0;
863:     dPD_dVr = dPD_dVi = dQD_dVr = dQD_dVi = 0.0;
864:     for (k=0; k < ld_nsegsp[i]; k++) {
865:       PD      += ld_alphap[k]*PD0[i]*PetscPowScalar((Vm/Vm0),ld_betap[k]);
866:       dPD_dVr += ld_alphap[k]*ld_betap[k]*PD0[i]*PetscPowScalar((1/Vm0),ld_betap[k])*Vr*PetscPowScalar(Vm,(ld_betap[k]-2));
867:       dPD_dVi += ld_alphap[k]*ld_betap[k]*PD0[i]*PetscPowScalar((1/Vm0),ld_betap[k])*Vi*PetscPowScalar(Vm,(ld_betap[k]-2));
868:     }
869:     for (k=0; k < ld_nsegsq[i]; k++) {
870:       QD      += ld_alphaq[k]*QD0[i]*PetscPowScalar((Vm/Vm0),ld_betaq[k]);
871:       dQD_dVr += ld_alphaq[k]*ld_betaq[k]*QD0[i]*PetscPowScalar((1/Vm0),ld_betaq[k])*Vr*PetscPowScalar(Vm,(ld_betaq[k]-2));
872:       dQD_dVi += ld_alphaq[k]*ld_betaq[k]*QD0[i]*PetscPowScalar((1/Vm0),ld_betaq[k])*Vi*PetscPowScalar(Vm,(ld_betaq[k]-2));
873:     }

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

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

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

884:     /*    fnet[2*lbus[i]]   += IDi; */
885:     row[0] = net_start + 2*lbus[i];
886:     col[0] = net_start + 2*lbus[i];  col[1] = net_start + 2*lbus[i]+1;
887:     val[0] = dIDi_dVr;               val[1] = dIDi_dVi;
888:     MatSetValues(J,1,row,2,col,val,ADD_VALUES);
889:     /*    fnet[2*lbus[i]+1] += IDr; */
890:     row[0] = net_start + 2*lbus[i]+1;
891:     col[0] = net_start + 2*lbus[i];  col[1] = net_start + 2*lbus[i]+1;
892:     val[0] = dIDr_dVr;               val[1] = dIDr_dVi;
893:     MatSetValues(J,1,row,2,col,val,ADD_VALUES);
894:   }
895:   VecRestoreArrayRead(user->V0,&v0);

897:   VecRestoreArrayRead(Xgen,&xgen);
898:   VecRestoreArrayRead(Xnet,&xnet);

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

902:   MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
903:   MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
904:   return 0;
905: }

907: /*
908:    J = [I, 0
909:         dg_dx, dg_dy]
910: */
911: PetscErrorCode AlgJacobian(SNES snes,Vec X,Mat A,Mat B,void *ctx)
912: {
913:   Userctx        *user=(Userctx*)ctx;

915:   ResidualJacobian(X,A,B,ctx);
916:   MatSetOption(A,MAT_KEEP_NONZERO_PATTERN,PETSC_TRUE);
917:   MatZeroRowsIS(A,user->is_diff,1.0,NULL,NULL);
918:   return 0;
919: }

921: /*
922:    J = [-df_dx, -df_dy
923:         dg_dx, dg_dy]
924: */

926: PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec X,Mat A,Mat B,void *ctx)
927: {
928:   Userctx        *user=(Userctx*)ctx;

930:   user->t = t;

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

934:   return 0;
935: }

937: /*
938:    J = [df_dx-aI, df_dy
939:         dg_dx, dg_dy]
940: */

942: PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,Userctx *user)
943: {
944:   PetscScalar    atmp = (PetscScalar) a;
945:   PetscInt       i,row;

947:   user->t = t;

949:   RHSJacobian(ts,t,X,A,B,user);
950:   MatScale(B,-1.0);
951:   for (i=0;i < ngen;i++) {
952:     row = 9*i;
953:     MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
954:     row  = 9*i+1;
955:     MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
956:     row  = 9*i+2;
957:     MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
958:     row  = 9*i+3;
959:     MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
960:     row  = 9*i+6;
961:     MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
962:     row  = 9*i+7;
963:     MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
964:     row  = 9*i+8;
965:     MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
966:   }
967:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
968:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
969:   return 0;
970: }

972: int main(int argc,char **argv)
973: {
974:   TS               ts;
975:   SNES              snes_alg;
976:   PetscErrorCode    ierr;
977:   PetscMPIInt       size;
978:   Userctx           user;
979:   PetscViewer       Xview,Ybusview,viewer;
980:   Vec               X,F_alg;
981:   Mat               J,A;
982:   PetscInt          i,idx,*idx2;
983:   Vec               Xdot;
984:   PetscScalar       *x,*mat,*amat;
985:   const PetscScalar *rmat;
986:   Vec               vatol;
987:   PetscInt          *direction;
988:   PetscBool         *terminate;
989:   const PetscInt    *idx3;
990:   PetscScalar       *vatoli;
991:   PetscInt          k;

993:   PetscInitialize(&argc,&argv,"petscoptions",help);
994:   MPI_Comm_size(PETSC_COMM_WORLD,&size);

997:   user.neqs_gen   = 9*ngen; /* # eqs. for generator subsystem */
998:   user.neqs_net   = 2*nbus; /* # eqs. for network subsystem   */
999:   user.neqs_pgrid = user.neqs_gen + user.neqs_net;

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

1003:   PetscMalloc1(7*ngen,&idx2);
1004:   for (i=0; i<ngen; i++) {
1005:     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;
1006:     idx2[7*i+4] = 9*i+6; idx2[7*i+5] = 9*i+7; idx2[7*i+6] = 9*i+8;
1007:   }
1008:   ISCreateGeneral(PETSC_COMM_WORLD,7*ngen,idx2,PETSC_COPY_VALUES,&user.is_diff);
1009:   ISComplement(user.is_diff,0,user.neqs_pgrid,&user.is_alg);
1010:   PetscFree(idx2);

1012:   /* Read initial voltage vector and Ybus */
1013:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"X.bin",FILE_MODE_READ,&Xview);
1014:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"Ybus.bin",FILE_MODE_READ,&Ybusview);

1016:   VecCreate(PETSC_COMM_WORLD,&user.V0);
1017:   VecSetSizes(user.V0,PETSC_DECIDE,user.neqs_net);
1018:   VecLoad(user.V0,Xview);

1020:   MatCreate(PETSC_COMM_WORLD,&user.Ybus);
1021:   MatSetSizes(user.Ybus,PETSC_DECIDE,PETSC_DECIDE,user.neqs_net,user.neqs_net);
1022:   MatSetType(user.Ybus,MATBAIJ);
1023:   /*  MatSetBlockSize(user.Ybus,2); */
1024:   MatLoad(user.Ybus,Ybusview);

1026:   /* Set run time options */
1027:   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Transient stability fault options","");
1028:   {
1029:     user.tfaulton  = 1.0;
1030:     user.tfaultoff = 1.2;
1031:     user.Rfault    = 0.0001;
1032:     user.setisdiff = PETSC_FALSE;
1033:     user.semiexplicit = PETSC_FALSE;
1034:     user.faultbus  = 8;
1035:     PetscOptionsReal("-tfaulton","","",user.tfaulton,&user.tfaulton,NULL);
1036:     PetscOptionsReal("-tfaultoff","","",user.tfaultoff,&user.tfaultoff,NULL);
1037:     PetscOptionsInt("-faultbus","","",user.faultbus,&user.faultbus,NULL);
1038:     user.t0        = 0.0;
1039:     user.tmax      = 5.0;
1040:     PetscOptionsReal("-t0","","",user.t0,&user.t0,NULL);
1041:     PetscOptionsReal("-tmax","","",user.tmax,&user.tmax,NULL);
1042:     PetscOptionsBool("-setisdiff","","",user.setisdiff,&user.setisdiff,NULL);
1043:     PetscOptionsBool("-dae_semiexplicit","","",user.semiexplicit,&user.semiexplicit,NULL);
1044:   }
1045:   PetscOptionsEnd();

1047:   PetscViewerDestroy(&Xview);
1048:   PetscViewerDestroy(&Ybusview);

1050:   /* Create DMs for generator and network subsystems */
1051:   DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,user.neqs_gen,1,1,NULL,&user.dmgen);
1052:   DMSetOptionsPrefix(user.dmgen,"dmgen_");
1053:   DMSetFromOptions(user.dmgen);
1054:   DMSetUp(user.dmgen);
1055:   DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,user.neqs_net,1,1,NULL,&user.dmnet);
1056:   DMSetOptionsPrefix(user.dmnet,"dmnet_");
1057:   DMSetFromOptions(user.dmnet);
1058:   DMSetUp(user.dmnet);
1059:   /* Create a composite DM packer and add the two DMs */
1060:   DMCompositeCreate(PETSC_COMM_WORLD,&user.dmpgrid);
1061:   DMSetOptionsPrefix(user.dmpgrid,"pgrid_");
1062:   DMCompositeAddDM(user.dmpgrid,user.dmgen);
1063:   DMCompositeAddDM(user.dmpgrid,user.dmnet);

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

1067:   MatCreate(PETSC_COMM_WORLD,&J);
1068:   MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,user.neqs_pgrid,user.neqs_pgrid);
1069:   MatSetFromOptions(J);
1070:   PreallocateJacobian(J,&user);

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

1075:   MatCreateSeqDense(PETSC_COMM_SELF,user.neqs_pgrid+1,1002,NULL,&user.Sol);
1076:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1077:      Create timestepping solver context
1078:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1079:   TSCreate(PETSC_COMM_WORLD,&ts);
1080:   TSSetProblemType(ts,TS_NONLINEAR);
1081:   if (user.semiexplicit) {
1082:     TSSetType(ts,TSRK);
1083:     TSSetRHSFunction(ts,NULL,RHSFunction,&user);
1084:     TSSetRHSJacobian(ts,J,J,RHSJacobian,&user);
1085:   } else {
1086:     TSSetType(ts,TSCN);
1087:     TSSetEquationType(ts,TS_EQ_DAE_IMPLICIT_INDEX1);
1088:     TSSetIFunction(ts,NULL,(TSIFunction) IFunction,&user);
1089:     TSSetIJacobian(ts,J,J,(TSIJacobian)IJacobian,&user);
1090:   }
1091:   TSSetApplicationContext(ts,&user);

1093:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1094:      Set initial conditions
1095:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1096:   SetInitialGuess(X,&user);
1097:   /* Just to set up the Jacobian structure */

1099:   VecDuplicate(X,&Xdot);
1100:   IJacobian(ts,0.0,X,Xdot,0.0,J,J,&user);
1101:   VecDestroy(&Xdot);

1103:   /* Save initial solution */

1105:   idx=user.stepnum*(user.neqs_pgrid+1);
1106:   MatDenseGetArray(user.Sol,&mat);
1107:   VecGetArray(X,&x);

1109:   mat[idx] = 0.0;

1111:   PetscArraycpy(mat+idx+1,x,user.neqs_pgrid);
1112:   MatDenseRestoreArray(user.Sol,&mat);
1113:   VecRestoreArray(X,&x);
1114:   user.stepnum++;

1116:   TSSetMaxTime(ts,user.tmax);
1117:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);
1118:   TSSetTimeStep(ts,0.01);
1119:   TSSetFromOptions(ts);
1120:   TSSetPostStep(ts,SaveSolution);
1121:   TSSetSolution(ts,X);

1123:   PetscMalloc1((2*ngen+2),&direction);
1124:   PetscMalloc1((2*ngen+2),&terminate);
1125:   direction[0] = direction[1] = 1;
1126:   terminate[0] = terminate[1] = PETSC_FALSE;
1127:   for (i=0; i < ngen;i++) {
1128:     direction[2+2*i] = -1; direction[2+2*i+1] = 1;
1129:     terminate[2+2*i] = terminate[2+2*i+1] = PETSC_FALSE;
1130:   }

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

1134:   if (user.semiexplicit) {
1135:     /* Use a semi-explicit approach with the time-stepping done by an explicit method and the
1136:        algrebraic part solved via PostStage and PostEvaluate callbacks
1137:     */
1138:     TSSetType(ts,TSRK);
1139:     TSSetPostStage(ts,PostStage);
1140:     TSSetPostEvaluate(ts,PostEvaluate);
1141:   }

1143:   if (user.setisdiff) {
1144:     /* Create vector of absolute tolerances and set the algebraic part to infinity */
1145:     VecDuplicate(X,&vatol);
1146:     VecSet(vatol,100000.0);
1147:     VecGetArray(vatol,&vatoli);
1148:     ISGetIndices(user.is_diff,&idx3);
1149:     for (k=0; k < 7*ngen; k++) vatoli[idx3[k]] = 1e-2;
1150:     VecRestoreArray(vatol,&vatoli);
1151:   }

1153:   /* Create the nonlinear solver for solving the algebraic system */
1154:   /* Note that although the algebraic system needs to be solved only for
1155:      Idq and V, we reuse the entire system including xgen. The xgen
1156:      variables are held constant by setting their residuals to 0 and
1157:      putting a 1 on the Jacobian diagonal for xgen rows
1158:   */

1160:   VecDuplicate(X,&F_alg);
1161:   SNESCreate(PETSC_COMM_WORLD,&snes_alg);
1162:   SNESSetFunction(snes_alg,F_alg,AlgFunction,&user);
1163:   SNESSetJacobian(snes_alg,J,J,AlgJacobian,&user);
1164:   SNESSetFromOptions(snes_alg);

1166:   user.snes_alg=snes_alg;

1168:   /* Solve */
1169:   TSSolve(ts,X);

1171:   MatAssemblyBegin(user.Sol,MAT_FINAL_ASSEMBLY);
1172:   MatAssemblyEnd(user.Sol,MAT_FINAL_ASSEMBLY);

1174:   MatCreateSeqDense(PETSC_COMM_SELF,user.neqs_pgrid+1,user.stepnum,NULL,&A);
1175:   MatDenseGetArrayRead(user.Sol,&rmat);
1176:   MatDenseGetArray(A,&amat);
1177:   PetscArraycpy(amat,rmat,user.stepnum*(user.neqs_pgrid+1));
1178:   MatDenseRestoreArray(A,&amat);
1179:   MatDenseRestoreArrayRead(user.Sol,&rmat);
1180:   PetscViewerBinaryOpen(PETSC_COMM_SELF,"out.bin",FILE_MODE_WRITE,&viewer);
1181:   MatView(A,viewer);
1182:   PetscViewerDestroy(&viewer);
1183:   MatDestroy(&A);

1185:   PetscFree(direction);
1186:   PetscFree(terminate);
1187:   SNESDestroy(&snes_alg);
1188:   VecDestroy(&F_alg);
1189:   MatDestroy(&J);
1190:   MatDestroy(&user.Ybus);
1191:   MatDestroy(&user.Sol);
1192:   VecDestroy(&X);
1193:   VecDestroy(&user.V0);
1194:   DMDestroy(&user.dmgen);
1195:   DMDestroy(&user.dmnet);
1196:   DMDestroy(&user.dmpgrid);
1197:   ISDestroy(&user.is_diff);
1198:   ISDestroy(&user.is_alg);
1199:   TSDestroy(&ts);
1200:   if (user.setisdiff) {
1201:     VecDestroy(&vatol);
1202:   }
1203:   PetscFinalize();
1204:   return 0;
1205: }

1207: /*TEST

1209:    build:
1210:       requires: double !complex !defined(PETSC_USE_64BIT_INDICES)

1212:    test:
1213:       suffix: implicit
1214:       args: -ts_monitor -snes_monitor_short
1215:       localrunfiles: petscoptions X.bin Ybus.bin

1217:    test:
1218:       suffix: semiexplicit
1219:       args: -ts_monitor -snes_monitor_short -dae_semiexplicit -ts_rk_type 2a
1220:       localrunfiles: petscoptions X.bin Ybus.bin

1222:    test:
1223:       suffix: steprestart
1224:       # needs ARKIMEX methods with all implicit stages since the mass matrix is not the identity
1225:       args: -ts_monitor -snes_monitor_short -ts_type arkimex -ts_arkimex_type prssp2
1226:       localrunfiles: petscoptions X.bin Ybus.bin

1228: TEST*/