Actual source code: ex9bus.c

petsc-3.7.7 2017-09-25
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) */

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

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

107: typedef struct {
108:   DM          dmgen, dmnet; /* DMs to manage generator and network subsystem */
109:   DM          dmpgrid; /* Composite DM to manage the entire power grid */
110:   Mat         Ybus; /* Network admittance matrix */
111:   Vec         V0;  /* Initial voltage vector (Power flow solution) */
112:   PetscReal   tfaulton,tfaultoff; /* Fault on and off times */
113:   PetscInt    faultbus; /* Fault bus */
114:   PetscScalar Rfault;
115:   PetscReal   t0,tmax;
116:   PetscInt    neqs_gen,neqs_net,neqs_pgrid;
117:   Mat         Sol; /* Matrix to save solution at each time step */
118:   PetscInt    stepnum;
119:   PetscBool   alg_flg;
120:   PetscReal   t;
121:   IS          is_diff; /* indices for differential equations */
122:   IS          is_alg; /* indices for algebraic equations */
123:   PetscBool   setisdiff; /* TS computes truncation error based only on the differential variables */
124: } Userctx;


127: /* Converts from machine frame (dq) to network (phase a real,imag) reference frame */
130: PetscErrorCode dq2ri(PetscScalar Fd,PetscScalar Fq,PetscScalar delta,PetscScalar *Fr, PetscScalar *Fi)
131: {
133:   *Fr =  Fd*PetscSinScalar(delta) + Fq*PetscCosScalar(delta);
134:   *Fi = -Fd*PetscCosScalar(delta) + Fq*PetscSinScalar(delta);
135:   return(0);
136: }

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

149: /* Saves the solution at each time to a matrix */
152: PetscErrorCode SaveSolution(TS ts)
153: {
154:   PetscErrorCode    ierr;
155:   Userctx           *user;
156:   Vec               X;
157:   const PetscScalar *x;
158:   PetscScalar       *mat;
159:   PetscInt          idx;
160:   PetscReal         t;

163:   TSGetApplicationContext(ts,&user);
164:   TSGetTime(ts,&t);
165:   TSGetSolution(ts,&X);
166:   idx      = user->stepnum*(user->neqs_pgrid+1);
167:   MatDenseGetArray(user->Sol,&mat);
168:   VecGetArrayRead(X,&x);
169:   mat[idx] = t;
170:   PetscMemcpy(mat+idx+1,x,user->neqs_pgrid*sizeof(PetscScalar));
171:   MatDenseRestoreArray(user->Sol,&mat);
172:   VecRestoreArrayRead(X,&x);
173:   user->stepnum++;
174:   return(0);
175: }

179: PetscErrorCode SetInitialGuess(Vec X,Userctx *user)
180: {
182:   Vec            Xgen,Xnet;
183:   PetscScalar    *xgen,*xnet;
184:   PetscInt       i,idx=0;
185:   PetscScalar    Vr,Vi,IGr,IGi,Vm,Vm2;
186:   PetscScalar    Eqp,Edp,delta;
187:   PetscScalar    Efd,RF,VR; /* Exciter variables */
188:   PetscScalar    Id,Iq;  /* Generator dq axis currents */
189:   PetscScalar    theta,Vd,Vq,SE;

192:   M[0] = 2*H[0]/w_s; M[1] = 2*H[1]/w_s; M[2] = 2*H[2]/w_s;
193:   D[0] = 0.1*M[0]; D[1] = 0.1*M[1]; D[2] = 0.1*M[2];

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

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

200:   /* Generator subsystem initialization */
201:   VecGetArray(Xgen,&xgen);
202:   VecGetArray(Xnet,&xnet);

204:   for (i=0; i < ngen; i++) {
205:     Vr  = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */
206:     Vi  = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */
207:     Vm  = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm;
208:     IGr = (Vr*PG[i] + Vi*QG[i])/Vm2;
209:     IGi = (Vi*PG[i] - Vr*QG[i])/Vm2;

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

213:     theta = PETSC_PI/2.0 - delta;

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

218:     Vd = Vr*PetscCosScalar(theta) - Vi*PetscSinScalar(theta);
219:     Vq = Vr*PetscSinScalar(theta) + Vi*PetscCosScalar(theta);

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

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

226:     /* The generator variables are ordered as [Eqp,Edp,delta,w,Id,Iq] */
227:     xgen[idx]   = Eqp;
228:     xgen[idx+1] = Edp;
229:     xgen[idx+2] = delta;
230:     xgen[idx+3] = w_s;

232:     idx = idx + 4;

234:     xgen[idx]   = Id;
235:     xgen[idx+1] = Iq;

237:     idx = idx + 2;

239:     /* Exciter */
240:     Efd = Eqp + (Xd[i] - Xdp[i])*Id;
241:     SE  = k1[i]*PetscExpScalar(k2[i]*Efd);
242:     VR  =  KE[i]*Efd + SE;
243:     RF  =  KF[i]*Efd/TF[i];

245:     xgen[idx]   = Efd;
246:     xgen[idx+1] = RF;
247:     xgen[idx+2] = VR;

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

251:     idx = idx + 3;
252:   }

254:   VecRestoreArray(Xgen,&xgen);
255:   VecRestoreArray(Xnet,&xnet);

257:   /* VecView(Xgen,0); */
258:   DMCompositeGather(user->dmpgrid,X,INSERT_VALUES,Xgen,Xnet);
259:   DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);
260:   return(0);
261: }

263: /* Computes F = [f(x,y);g(x,y)] */
266: PetscErrorCode ResidualFunction(SNES snes,Vec X, Vec F, Userctx *user)
267: {
269:   Vec            Xgen,Xnet,Fgen,Fnet;
270:   PetscScalar    *xgen,*xnet,*fgen,*fnet;
271:   PetscInt       i,idx=0;
272:   PetscScalar    Vr,Vi,Vm,Vm2;
273:   PetscScalar    Eqp,Edp,delta,w; /* Generator variables */
274:   PetscScalar    Efd,RF,VR; /* Exciter variables */
275:   PetscScalar    Id,Iq;  /* Generator dq axis currents */
276:   PetscScalar    Vd,Vq,SE;
277:   PetscScalar    IGr,IGi,IDr,IDi;
278:   PetscScalar    Zdq_inv[4],det;
279:   PetscScalar    PD,QD,Vm0,*v0;
280:   PetscInt       k;

283:   VecZeroEntries(F);
284:   DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
285:   DMCompositeGetLocalVectors(user->dmpgrid,&Fgen,&Fnet);
286:   DMCompositeScatter(user->dmpgrid,X,Xgen,Xnet);
287:   DMCompositeScatter(user->dmpgrid,F,Fgen,Fnet);

289:   /* Network current balance residual IG + Y*V + IL = 0. Only YV is added here.
290:      The generator current injection, IG, and load current injection, ID are added later
291:   */
292:   /* Note that the values in Ybus are stored assuming the imaginary current balance
293:      equation is ordered first followed by real current balance equation for each bus.
294:      Thus imaginary current contribution goes in location 2*i, and
295:      real current contribution in 2*i+1
296:   */
297:   MatMult(user->Ybus,Xnet,Fnet);

299:   VecGetArray(Xgen,&xgen);
300:   VecGetArray(Xnet,&xnet);
301:   VecGetArray(Fgen,&fgen);
302:   VecGetArray(Fnet,&fnet);

304:   /* Generator subsystem */
305:   for (i=0; i < ngen; i++) {
306:     Eqp   = xgen[idx];
307:     Edp   = xgen[idx+1];
308:     delta = xgen[idx+2];
309:     w     = xgen[idx+3];
310:     Id    = xgen[idx+4];
311:     Iq    = xgen[idx+5];
312:     Efd   = xgen[idx+6];
313:     RF    = xgen[idx+7];
314:     VR    = xgen[idx+8];

316:     /* Generator differential equations */
317:     fgen[idx]   = (Eqp + (Xd[i] - Xdp[i])*Id - Efd)/Td0p[i];
318:     fgen[idx+1] = (Edp - (Xq[i] - Xqp[i])*Iq)/Tq0p[i];
319:     fgen[idx+2] = -w + w_s;
320:     fgen[idx+3] = (-TM[i] + Edp*Id + Eqp*Iq + (Xqp[i] - Xdp[i])*Id*Iq + D[i]*(w - w_s))/M[i];

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

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

329:     Zdq_inv[0] = Rs[i]/det;
330:     Zdq_inv[1] = Xqp[i]/det;
331:     Zdq_inv[2] = -Xdp[i]/det;
332:     Zdq_inv[3] = Rs[i]/det;

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

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

340:     fnet[2*gbus[i]]   -= IGi;
341:     fnet[2*gbus[i]+1] -= IGr;

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

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

347:     /* Exciter differential equations */
348:     fgen[idx+6] = (KE[i]*Efd + SE - VR)/TE[i];
349:     fgen[idx+7] = (RF - KF[i]*Efd/TF[i])/TF[i];
350:     fgen[idx+8] = (VR - KA[i]*RF + KA[i]*KF[i]*Efd/TF[i] - KA[i]*(Vref[i] - Vm))/TA[i];

352:     idx = idx + 9;
353:   }

355:   VecGetArray(user->V0,&v0);
356:   for (i=0; i < nload; i++) {
357:     Vr  = xnet[2*lbus[i]]; /* Real part of load bus voltage */
358:     Vi  = xnet[2*lbus[i]+1]; /* Imaginary part of the load bus voltage */
359:     Vm  = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm;
360:     Vm0 = PetscSqrtScalar(v0[2*lbus[i]]*v0[2*lbus[i]] + v0[2*lbus[i]+1]*v0[2*lbus[i]+1]);
361:     PD  = QD = 0.0;
362:     for (k=0; k < ld_nsegsp[i]; k++) PD += ld_alphap[k]*PD0[i]*PetscPowScalar((Vm/Vm0),ld_betap[k]);
363:     for (k=0; k < ld_nsegsq[i]; k++) QD += ld_alphaq[k]*QD0[i]*PetscPowScalar((Vm/Vm0),ld_betaq[k]);

365:     /* Load currents */
366:     IDr = (PD*Vr + QD*Vi)/Vm2;
367:     IDi = (-QD*Vr + PD*Vi)/Vm2;

369:     fnet[2*lbus[i]]   += IDi;
370:     fnet[2*lbus[i]+1] += IDr;
371:   }
372:   VecRestoreArray(user->V0,&v0);

374:   VecRestoreArray(Xgen,&xgen);
375:   VecRestoreArray(Xnet,&xnet);
376:   VecRestoreArray(Fgen,&fgen);
377:   VecRestoreArray(Fnet,&fnet);

379:   DMCompositeGather(user->dmpgrid,F,INSERT_VALUES,Fgen,Fnet);
380:   DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);
381:   DMCompositeRestoreLocalVectors(user->dmpgrid,&Fgen,&Fnet);
382:   return(0);
383: }

385: /* \dot{x} - f(x,y)
386:      g(x,y) = 0
387:  */
390: PetscErrorCode IFunction(TS ts,PetscReal t, Vec X, Vec Xdot, Vec F, Userctx *user)
391: {
393:   SNES           snes;
394:   PetscScalar    *f,*xdot;
395:   PetscInt       i;

398:   user->t = t;

400:   TSGetSNES(ts,&snes);
401:   ResidualFunction(snes,X,F,user);
402:   VecGetArray(F,&f);
403:   VecGetArray(Xdot,&xdot);
404:   for (i=0;i < ngen;i++) {
405:     f[9*i]   += xdot[9*i];
406:     f[9*i+1] += xdot[9*i+1];
407:     f[9*i+2] += xdot[9*i+2];
408:     f[9*i+3] += xdot[9*i+3];
409:     f[9*i+6] += xdot[9*i+6];
410:     f[9*i+7] += xdot[9*i+7];
411:     f[9*i+8] += xdot[9*i+8];
412:   }
413:   VecRestoreArray(F,&f);
414:   VecRestoreArray(Xdot,&xdot);
415:   return(0);
416: }

418: /* This function is used for solving the algebraic system only during fault on and
419:    off times. It computes the entire F and then zeros out the part corresponding to
420:    differential equations
421:  F = [0;g(y)];
422: */
425: PetscErrorCode AlgFunction(SNES snes, Vec X, Vec F, void *ctx)
426: {
428:   Userctx        *user=(Userctx*)ctx;
429:   PetscScalar    *f;
430:   PetscInt       i;

433:   ResidualFunction(snes,X,F,user);
434:   VecGetArray(F,&f);
435:   for (i=0; i < ngen; i++) {
436:     f[9*i]   = 0;
437:     f[9*i+1] = 0;
438:     f[9*i+2] = 0;
439:     f[9*i+3] = 0;
440:     f[9*i+6] = 0;
441:     f[9*i+7] = 0;
442:     f[9*i+8] = 0;
443:   }
444:   VecRestoreArray(F,&f);
445:   return(0);
446: }

450: PetscErrorCode PreallocateJacobian(Mat J, Userctx *user)
451: {
453:   PetscInt       *d_nnz;
454:   PetscInt       i,idx=0,start=0;
455:   PetscInt       ncols;

458:   PetscMalloc1(user->neqs_pgrid,&d_nnz);
459:   for (i=0; i<user->neqs_pgrid; i++) d_nnz[i] = 0;
460:   /* Generator subsystem */
461:   for (i=0; i < ngen; i++) {

463:     d_nnz[idx]   += 3;
464:     d_nnz[idx+1] += 2;
465:     d_nnz[idx+2] += 2;
466:     d_nnz[idx+3] += 5;
467:     d_nnz[idx+4] += 6;
468:     d_nnz[idx+5] += 6;

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

473:     d_nnz[idx+6] += 2;
474:     d_nnz[idx+7] += 2;
475:     d_nnz[idx+8] += 5;

477:     idx = idx + 9;
478:   }

480:   start = user->neqs_gen;

482:   for (i=0; i < nbus; i++) {
483:     MatGetRow(user->Ybus,2*i,&ncols,NULL,NULL);
484:     d_nnz[start+2*i]   += ncols;
485:     d_nnz[start+2*i+1] += ncols;
486:     MatRestoreRow(user->Ybus,2*i,&ncols,NULL,NULL);
487:   }

489:   MatSeqAIJSetPreallocation(J,0,d_nnz);

491:   PetscFree(d_nnz);
492:   return(0);
493: }

495: /*
496:    J = [-df_dx, -df_dy
497:         dg_dx, dg_dy]
498: */
501: PetscErrorCode ResidualJacobian(SNES snes,Vec X,Mat J,Mat B,void *ctx)
502: {
504:   Userctx        *user=(Userctx*)ctx;
505:   Vec            Xgen,Xnet;
506:   PetscScalar    *xgen,*xnet;
507:   PetscInt       i,idx=0;
508:   PetscScalar    Vr,Vi,Vm,Vm2;
509:   PetscScalar    Eqp,Edp,delta; /* Generator variables */
510:   PetscScalar    Efd;
511:   PetscScalar    Id,Iq;  /* Generator dq axis currents */
512:   PetscScalar    Vd,Vq;
513:   PetscScalar    val[10];
514:   PetscInt       row[2],col[10];
515:   PetscInt       net_start=user->neqs_gen;
516:   PetscScalar    Zdq_inv[4],det;
517:   PetscScalar    dVd_dVr,dVd_dVi,dVq_dVr,dVq_dVi,dVd_ddelta,dVq_ddelta;
518:   PetscScalar    dIGr_ddelta,dIGi_ddelta,dIGr_dId,dIGr_dIq,dIGi_dId,dIGi_dIq;
519:   PetscScalar    dSE_dEfd;
520:   PetscScalar    dVm_dVd,dVm_dVq,dVm_dVr,dVm_dVi;
521:   PetscInt          ncols;
522:   const PetscInt    *cols;
523:   const PetscScalar *yvals;
524:   PetscInt          k;
525:   PetscScalar PD,QD,Vm0,*v0,Vm4;
526:   PetscScalar dPD_dVr,dPD_dVi,dQD_dVr,dQD_dVi;
527:   PetscScalar dIDr_dVr,dIDr_dVi,dIDi_dVr,dIDi_dVi;


531:   MatZeroEntries(B);
532:   DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
533:   DMCompositeScatter(user->dmpgrid,X,Xgen,Xnet);

535:   VecGetArray(Xgen,&xgen);
536:   VecGetArray(Xnet,&xnet);

538:   /* Generator subsystem */
539:   for (i=0; i < ngen; i++) {
540:     Eqp   = xgen[idx];
541:     Edp   = xgen[idx+1];
542:     delta = xgen[idx+2];
543:     Id    = xgen[idx+4];
544:     Iq    = xgen[idx+5];
545:     Efd   = xgen[idx+6];

547:     /*    fgen[idx]   = (Eqp + (Xd[i] - Xdp[i])*Id - Efd)/Td0p[i]; */
548:     row[0] = idx;
549:     col[0] = idx;           col[1] = idx+4;          col[2] = idx+6;
550:     val[0] = 1/ Td0p[i]; val[1] = (Xd[i] - Xdp[i])/ Td0p[i]; val[2] = -1/Td0p[i];

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

554:     /*    fgen[idx+1] = (Edp - (Xq[i] - Xqp[i])*Iq)/Tq0p[i]; */
555:     row[0] = idx + 1;
556:     col[0] = idx + 1;       col[1] = idx+5;
557:     val[0] = 1/Tq0p[i]; val[1] = -(Xq[i] - Xqp[i])/Tq0p[i];
558:     MatSetValues(J,1,row,2,col,val,INSERT_VALUES);

560:     /*    fgen[idx+2] = - w + w_s; */
561:     row[0] = idx + 2;
562:     col[0] = idx + 2; col[1] = idx + 3;
563:     val[0] = 0;       val[1] = -1;
564:     MatSetValues(J,1,row,2,col,val,INSERT_VALUES);

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

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

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

578:     Zdq_inv[0] = Rs[i]/det;
579:     Zdq_inv[1] = Xqp[i]/det;
580:     Zdq_inv[2] = -Xdp[i]/det;
581:     Zdq_inv[3] = Rs[i]/det;

583:     dVd_dVr    = PetscSinScalar(delta); dVd_dVi = -PetscCosScalar(delta);
584:     dVq_dVr    = PetscCosScalar(delta); dVq_dVi = PetscSinScalar(delta);
585:     dVd_ddelta = Vr*PetscCosScalar(delta) + Vi*PetscSinScalar(delta);
586:     dVq_ddelta = -Vr*PetscSinScalar(delta) + Vi*PetscCosScalar(delta);

588:     /*    fgen[idx+4] = Zdq_inv[0]*(-Edp + Vd) + Zdq_inv[1]*(-Eqp + Vq) + Id; */
589:     row[0] = idx+4;
590:     col[0] = idx;         col[1] = idx+1;        col[2] = idx + 2;
591:     val[0] = -Zdq_inv[1]; val[1] = -Zdq_inv[0];  val[2] = Zdq_inv[0]*dVd_ddelta + Zdq_inv[1]*dVq_ddelta;
592:     col[3] = idx + 4; col[4] = net_start+2*gbus[i];                     col[5] = net_start + 2*gbus[i]+1;
593:     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;
594:     MatSetValues(J,1,row,6,col,val,INSERT_VALUES);

596:     /*  fgen[idx+5] = Zdq_inv[2]*(-Edp + Vd) + Zdq_inv[3]*(-Eqp + Vq) + Iq; */
597:     row[0] = idx+5;
598:     col[0] = idx;         col[1] = idx+1;        col[2] = idx + 2;
599:     val[0] = -Zdq_inv[3]; val[1] = -Zdq_inv[2];  val[2] = Zdq_inv[2]*dVd_ddelta + Zdq_inv[3]*dVq_ddelta;
600:     col[3] = idx + 5; col[4] = net_start+2*gbus[i];                     col[5] = net_start + 2*gbus[i]+1;
601:     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;
602:     MatSetValues(J,1,row,6,col,val,INSERT_VALUES);

604:     dIGr_ddelta = Id*PetscCosScalar(delta) - Iq*PetscSinScalar(delta);
605:     dIGi_ddelta = Id*PetscSinScalar(delta) + Iq*PetscCosScalar(delta);
606:     dIGr_dId    = PetscSinScalar(delta);  dIGr_dIq = PetscCosScalar(delta);
607:     dIGi_dId    = -PetscCosScalar(delta); dIGi_dIq = PetscSinScalar(delta);

609:     /* fnet[2*gbus[i]]   -= IGi; */
610:     row[0] = net_start + 2*gbus[i];
611:     col[0] = idx+2;        col[1] = idx + 4;   col[2] = idx + 5;
612:     val[0] = -dIGi_ddelta; val[1] = -dIGi_dId; val[2] = -dIGi_dIq;
613:     MatSetValues(J,1,row,3,col,val,INSERT_VALUES);

615:     /* fnet[2*gbus[i]+1]   -= IGr; */
616:     row[0] = net_start + 2*gbus[i]+1;
617:     col[0] = idx+2;        col[1] = idx + 4;   col[2] = idx + 5;
618:     val[0] = -dIGr_ddelta; val[1] = -dIGr_dId; val[2] = -dIGr_dIq;
619:     MatSetValues(J,1,row,3,col,val,INSERT_VALUES);

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

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

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

628:     row[0] = idx + 6;
629:     col[0] = idx + 6;                     col[1] = idx + 8;
630:     val[0] = (KE[i] + dSE_dEfd)/TE[i];  val[1] = -1/TE[i];
631:     MatSetValues(J,1,row,2,col,val,INSERT_VALUES);

633:     /* Exciter differential equations */

635:     /*    fgen[idx+7] = (RF - KF[i]*Efd/TF[i])/TF[i]; */
636:     row[0] = idx + 7;
637:     col[0] = idx + 6;       col[1] = idx + 7;
638:     val[0] = (-KF[i]/TF[i])/TF[i];  val[1] = 1/TF[i];
639:     MatSetValues(J,1,row,2,col,val,INSERT_VALUES);

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

644:     dVm_dVd    = Vd/Vm; dVm_dVq = Vq/Vm;
645:     dVm_dVr    = dVm_dVd*dVd_dVr + dVm_dVq*dVq_dVr;
646:     dVm_dVi    = dVm_dVd*dVd_dVi + dVm_dVq*dVq_dVi;
647:     row[0]     = idx + 8;
648:     col[0]     = idx + 6;           col[1] = idx + 7; col[2] = idx + 8;
649:     val[0]     = (KA[i]*KF[i]/TF[i])/TA[i]; val[1] = -KA[i]/TA[i];  val[2] = 1/TA[i];
650:     col[3]     = net_start + 2*gbus[i]; col[4] = net_start + 2*gbus[i]+1;
651:     val[3]     = KA[i]*dVm_dVr/TA[i];         val[4] = KA[i]*dVm_dVi/TA[i];
652:     MatSetValues(J,1,row,5,col,val,INSERT_VALUES);
653:     idx        = idx + 9;
654:   }

656:   for (i=0; i<nbus; i++) {
657:     MatGetRow(user->Ybus,2*i,&ncols,&cols,&yvals);
658:     row[0] = net_start + 2*i;
659:     for (k=0; k<ncols; k++) {
660:       col[k] = net_start + cols[k];
661:       val[k] = yvals[k];
662:     }
663:     MatSetValues(J,1,row,ncols,col,val,INSERT_VALUES);
664:     MatRestoreRow(user->Ybus,2*i,&ncols,&cols,&yvals);

666:     MatGetRow(user->Ybus,2*i+1,&ncols,&cols,&yvals);
667:     row[0] = net_start + 2*i+1;
668:     for (k=0; k<ncols; k++) {
669:       col[k] = net_start + cols[k];
670:       val[k] = yvals[k];
671:     }
672:     MatSetValues(J,1,row,ncols,col,val,INSERT_VALUES);
673:     MatRestoreRow(user->Ybus,2*i+1,&ncols,&cols,&yvals);
674:   }

676:   MatAssemblyBegin(J,MAT_FLUSH_ASSEMBLY);
677:   MatAssemblyEnd(J,MAT_FLUSH_ASSEMBLY);

679:   VecGetArray(user->V0,&v0);
680:   for (i=0; i < nload; i++) {
681:     Vr      = xnet[2*lbus[i]]; /* Real part of load bus voltage */
682:     Vi      = xnet[2*lbus[i]+1]; /* Imaginary part of the load bus voltage */
683:     Vm      = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm; Vm4 = Vm2*Vm2;
684:     Vm0     = PetscSqrtScalar(v0[2*lbus[i]]*v0[2*lbus[i]] + v0[2*lbus[i]+1]*v0[2*lbus[i]+1]);
685:     PD      = QD = 0.0;
686:     dPD_dVr = dPD_dVi = dQD_dVr = dQD_dVi = 0.0;
687:     for (k=0; k < ld_nsegsp[i]; k++) {
688:       PD      += ld_alphap[k]*PD0[i]*PetscPowScalar((Vm/Vm0),ld_betap[k]);
689:       dPD_dVr += ld_alphap[k]*ld_betap[k]*PD0[i]*PetscPowScalar((1/Vm0),ld_betap[k])*Vr*PetscPowScalar(Vm,(ld_betap[k]-2));
690:       dPD_dVi += ld_alphap[k]*ld_betap[k]*PD0[i]*PetscPowScalar((1/Vm0),ld_betap[k])*Vi*PetscPowScalar(Vm,(ld_betap[k]-2));
691:     }
692:     for (k=0; k < ld_nsegsq[i]; k++) {
693:       QD      += ld_alphaq[k]*QD0[i]*PetscPowScalar((Vm/Vm0),ld_betaq[k]);
694:       dQD_dVr += ld_alphaq[k]*ld_betaq[k]*QD0[i]*PetscPowScalar((1/Vm0),ld_betaq[k])*Vr*PetscPowScalar(Vm,(ld_betaq[k]-2));
695:       dQD_dVi += ld_alphaq[k]*ld_betaq[k]*QD0[i]*PetscPowScalar((1/Vm0),ld_betaq[k])*Vi*PetscPowScalar(Vm,(ld_betaq[k]-2));
696:     }

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

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

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


708:     /*    fnet[2*lbus[i]]   += IDi; */
709:     row[0] = net_start + 2*lbus[i];
710:     col[0] = net_start + 2*lbus[i];  col[1] = net_start + 2*lbus[i]+1;
711:     val[0] = dIDi_dVr;               val[1] = dIDi_dVi;
712:     MatSetValues(J,1,row,2,col,val,ADD_VALUES);
713:     /*    fnet[2*lbus[i]+1] += IDr; */
714:     row[0] = net_start + 2*lbus[i]+1;
715:     col[0] = net_start + 2*lbus[i];  col[1] = net_start + 2*lbus[i]+1;
716:     val[0] = dIDr_dVr;               val[1] = dIDr_dVi;
717:     MatSetValues(J,1,row,2,col,val,ADD_VALUES);
718:   }
719:   VecRestoreArray(user->V0,&v0);

721:   VecRestoreArray(Xgen,&xgen);
722:   VecRestoreArray(Xnet,&xnet);

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

726:   MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
727:   MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
728:   return(0);
729: }

731: /*
732:    J = [I, 0
733:         dg_dx, dg_dy]
734: */
737: PetscErrorCode AlgJacobian(SNES snes,Vec X,Mat A,Mat B,void *ctx)
738: {
740:   Userctx        *user=(Userctx*)ctx;

743:   ResidualJacobian(snes,X,A,B,ctx);
744:   MatSetOption(A,MAT_KEEP_NONZERO_PATTERN,PETSC_TRUE);
745:   MatZeroRowsIS(A,user->is_diff,1.0,NULL,NULL);
746:   return(0);
747: }

749: /*
750:    J = [a*I-df_dx, -df_dy
751:         dg_dx, dg_dy]
752: */

756: PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,Userctx *user)
757: {
759:   SNES           snes;
760:   PetscScalar    atmp = (PetscScalar) a;
761:   PetscInt       i,row;

764:   user->t = t;

766:   TSGetSNES(ts,&snes);
767:   ResidualJacobian(snes,X,A,B,user);
768:   for (i=0;i < ngen;i++) {
769:     row = 9*i;
770:     MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
771:     row  = 9*i+1;
772:     MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
773:     row  = 9*i+2;
774:     MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
775:     row  = 9*i+3;
776:     MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
777:     row  = 9*i+6;
778:     MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
779:     row  = 9*i+7;
780:     MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
781:     row  = 9*i+8;
782:     MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
783:   }
784:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
785:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
786:   return(0);
787: }

791: int main(int argc,char **argv)
792: {
793:   TS             ts;
794:   SNES           snes_alg;
796:   PetscMPIInt    size;
797:   Userctx        user;
798:   PetscViewer    Xview,Ybusview,viewer;
799:   Vec            X,F_alg;
800:   Mat            J,A;
801:   PetscInt       i,idx,*idx2,row_loc,col_loc;
802:   Vec            Xdot;
803:   PetscScalar    *x,*mat,val,*amat;
804:   Vec            vatol;

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

810:   user.neqs_gen   = 9*ngen; /* # eqs. for generator subsystem */
811:   user.neqs_net   = 2*nbus; /* # eqs. for network subsystem   */
812:   user.neqs_pgrid = user.neqs_gen + user.neqs_net;

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

816:   PetscMalloc1(7*ngen,&idx2);
817:   for (i=0; i<ngen; i++) {
818:     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;
819:     idx2[7*i+4] = 9*i+6; idx2[7*i+5] = 9*i+7; idx2[7*i+6] = 9*i+8;
820:   }
821:   ISCreateGeneral(PETSC_COMM_WORLD,7*ngen,idx2,PETSC_COPY_VALUES,&user.is_diff);
822:   ISComplement(user.is_diff,0,user.neqs_pgrid,&user.is_alg);
823:   PetscFree(idx2);

825:   /* Read initial voltage vector and Ybus */
826:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"X.bin",FILE_MODE_READ,&Xview);
827:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"Ybus.bin",FILE_MODE_READ,&Ybusview);

829:   VecCreate(PETSC_COMM_WORLD,&user.V0);
830:   VecSetSizes(user.V0,PETSC_DECIDE,user.neqs_net);
831:   VecLoad(user.V0,Xview);

833:   MatCreate(PETSC_COMM_WORLD,&user.Ybus);
834:   MatSetSizes(user.Ybus,PETSC_DECIDE,PETSC_DECIDE,user.neqs_net,user.neqs_net);
835:   MatSetType(user.Ybus,MATBAIJ);
836:   /*  MatSetBlockSize(user.Ybus,2); */
837:   MatLoad(user.Ybus,Ybusview);

839:   /* Set run time options */
840:   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Transient stability fault options","");
841:   {
842:     user.tfaulton  = 1.0;
843:     user.tfaultoff = 1.2;
844:     user.Rfault    = 0.0001;
845:     user.setisdiff = PETSC_FALSE;
846:     user.faultbus  = 8;
847:     PetscOptionsReal("-tfaulton","","",user.tfaulton,&user.tfaulton,NULL);
848:     PetscOptionsReal("-tfaultoff","","",user.tfaultoff,&user.tfaultoff,NULL);
849:     PetscOptionsInt("-faultbus","","",user.faultbus,&user.faultbus,NULL);
850:     user.t0        = 0.0;
851:     user.tmax      = 5.0;
852:     PetscOptionsReal("-t0","","",user.t0,&user.t0,NULL);
853:     PetscOptionsReal("-tmax","","",user.tmax,&user.tmax,NULL);
854:     PetscOptionsBool("-setisdiff","","",user.setisdiff,&user.setisdiff,NULL);
855:   }
856:   PetscOptionsEnd();

858:   PetscViewerDestroy(&Xview);
859:   PetscViewerDestroy(&Ybusview);

861:   /* Create DMs for generator and network subsystems */
862:   DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,user.neqs_gen,1,1,NULL,&user.dmgen);
863:   DMSetOptionsPrefix(user.dmgen,"dmgen_");
864:   DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,user.neqs_net,1,1,NULL,&user.dmnet);
865:   DMSetOptionsPrefix(user.dmnet,"dmnet_");
866:   /* Create a composite DM packer and add the two DMs */
867:   DMCompositeCreate(PETSC_COMM_WORLD,&user.dmpgrid);
868:   DMSetOptionsPrefix(user.dmpgrid,"pgrid_");
869:   DMCompositeAddDM(user.dmpgrid,user.dmgen);
870:   DMCompositeAddDM(user.dmpgrid,user.dmnet);

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

874:   MatCreate(PETSC_COMM_WORLD,&J);
875:   MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,user.neqs_pgrid,user.neqs_pgrid);
876:   MatSetFromOptions(J);
877:   PreallocateJacobian(J,&user);

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

882:   MatCreateSeqDense(PETSC_COMM_SELF,user.neqs_pgrid+1,1002,NULL,&user.Sol);
883:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
884:      Create timestepping solver context
885:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
886:   TSCreate(PETSC_COMM_WORLD,&ts);
887:   TSSetProblemType(ts,TS_NONLINEAR);
888:   TSSetEquationType(ts,TS_EQ_DAE_IMPLICIT_INDEX1);
889:   TSARKIMEXSetFullyImplicit(ts,PETSC_TRUE);
890:   TSSetIFunction(ts,NULL,(TSIFunction) IFunction,&user);
891:   TSSetIJacobian(ts,J,J,(TSIJacobian)IJacobian,&user);
892:   TSSetApplicationContext(ts,&user);

894:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
895:      Set initial conditions
896:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
897:   SetInitialGuess(X,&user);
898:   /* Just to set up the Jacobian structure */

900:   VecDuplicate(X,&Xdot);
901:   IJacobian(ts,0.0,X,Xdot,0.0,J,J,&user);
902:   VecDestroy(&Xdot);

904:   /* Save initial solution */

906:   idx=user.stepnum*(user.neqs_pgrid+1);
907:   MatDenseGetArray(user.Sol,&mat);
908:   VecGetArray(X,&x);

910:   mat[idx] = 0.0;

912:   PetscMemcpy(mat+idx+1,x,user.neqs_pgrid*sizeof(PetscScalar));
913:   MatDenseRestoreArray(user.Sol,&mat);
914:   VecRestoreArray(X,&x);
915:   user.stepnum++;

917:   TSSetDuration(ts,1000,user.tfaulton);
918:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
919:   TSSetInitialTimeStep(ts,0.0,0.01);
920:   TSSetFromOptions(ts);
921:   TSSetPostStep(ts,SaveSolution);

923:   if(user.setisdiff) {
924:     const PetscInt *idx;
925:     PetscScalar *vatoli;
926:     PetscInt k;
927:     /* Create vector of absolute tolerances and set the algebraic part to infinity */
928:     VecDuplicate(X,&vatol);
929:     VecSet(X,100000.0);
930:     VecGetArray(vatol,&vatoli);
931:     ISGetIndices(user.is_diff,&idx);
932:     for(k=0; k < 7*ngen; k++) vatoli[idx[k]] = 1e-2;
933:     VecRestoreArray(vatol,&vatoli);
934:   }

936:   user.alg_flg = PETSC_FALSE;
937:   /* Prefault period */
938:   TSSolve(ts,X);

940:   /* Create the nonlinear solver for solving the algebraic system */
941:   /* Note that although the algebraic system needs to be solved only for
942:      Idq and V, we reuse the entire system including xgen. The xgen
943:      variables are held constant by setting their residuals to 0 and
944:      putting a 1 on the Jacobian diagonal for xgen rows
945:   */

947:   VecDuplicate(X,&F_alg);
948:   SNESCreate(PETSC_COMM_WORLD,&snes_alg);
949:   SNESSetFunction(snes_alg,F_alg,AlgFunction,&user);
950:   MatZeroEntries(J);
951:   SNESSetJacobian(snes_alg,J,J,AlgJacobian,&user);
952:   SNESSetOptionsPrefix(snes_alg,"alg_");
953:   SNESSetFromOptions(snes_alg);

955:   /* Apply disturbance - resistive fault at user.faultbus */
956:   /* This is done by adding shunt conductance to the diagonal location
957:      in the Ybus matrix */
958:   row_loc = 2*user.faultbus; col_loc = 2*user.faultbus+1; /* Location for G */
959:   val     = 1/user.Rfault;
960:   MatSetValues(user.Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);
961:   row_loc = 2*user.faultbus+1; col_loc = 2*user.faultbus; /* Location for G */
962:   val     = 1/user.Rfault;
963:   MatSetValues(user.Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);

965:   MatAssemblyBegin(user.Ybus,MAT_FINAL_ASSEMBLY);
966:   MatAssemblyEnd(user.Ybus,MAT_FINAL_ASSEMBLY);

968:   user.alg_flg = PETSC_TRUE;
969:   /* Solve the algebraic equations */
970:   SNESSolve(snes_alg,NULL,X);

972:   /* Save fault-on solution */
973:   idx      = user.stepnum*(user.neqs_pgrid+1);
974:   MatDenseGetArray(user.Sol,&mat);
975:   VecGetArray(X,&x);
976:   mat[idx] = user.tfaulton;
977:   PetscMemcpy(mat+idx+1,x,user.neqs_pgrid*sizeof(PetscScalar));
978:   MatDenseRestoreArray(user.Sol,&mat);
979:   VecRestoreArray(X,&x);
980:   user.stepnum++;

982:   /* Disturbance period */
983:   TSSetDuration(ts,1000,user.tfaultoff);
984:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
985:   TSSetInitialTimeStep(ts,user.tfaulton,.01);

987:   user.alg_flg = PETSC_FALSE;

989:   TSSolve(ts,X);

991:   /* Remove the fault */
992:   row_loc = 2*user.faultbus; col_loc = 2*user.faultbus+1;
993:   val     = -1/user.Rfault;
994:   MatSetValues(user.Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);
995:   row_loc = 2*user.faultbus+1; col_loc = 2*user.faultbus;
996:   val     = -1/user.Rfault;
997:   MatSetValues(user.Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);

999:   MatAssemblyBegin(user.Ybus,MAT_FINAL_ASSEMBLY);
1000:   MatAssemblyEnd(user.Ybus,MAT_FINAL_ASSEMBLY);

1002:   MatZeroEntries(J);

1004:   user.alg_flg = PETSC_TRUE;

1006:   /* Solve the algebraic equations */
1007:   SNESSolve(snes_alg,NULL,X);

1009:   /* Save tfault off solution */
1010:   idx      = user.stepnum*(user.neqs_pgrid+1);
1011:   MatDenseGetArray(user.Sol,&mat);
1012:   VecGetArray(X,&x);
1013:   mat[idx] = user.tfaultoff;
1014:   PetscMemcpy(mat+idx+1,x,user.neqs_pgrid*sizeof(PetscScalar));
1015:   MatDenseRestoreArray(user.Sol,&mat);
1016:   VecRestoreArray(X,&x);
1017:   user.stepnum++;

1019:   /* Post-disturbance period */
1020:   TSSetDuration(ts,1000,user.tmax);
1021:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
1022:   TSSetInitialTimeStep(ts,user.tfaultoff,.01);

1024:   user.alg_flg = PETSC_TRUE;

1026:   TSSolve(ts,X);

1028:   MatAssemblyBegin(user.Sol,MAT_FINAL_ASSEMBLY);
1029:   MatAssemblyEnd(user.Sol,MAT_FINAL_ASSEMBLY);

1031:   MatCreateSeqDense(PETSC_COMM_SELF,user.neqs_pgrid+1,user.stepnum,NULL,&A);
1032:   MatDenseGetArray(user.Sol,&mat);
1033:   MatDenseGetArray(A,&amat);
1034:   PetscMemcpy(amat,mat,(user.stepnum*(user.neqs_pgrid+1))*sizeof(PetscScalar));
1035:   MatDenseRestoreArray(A,&amat);
1036:   MatDenseRestoreArray(user.Sol,&mat);
1037:   PetscViewerBinaryOpen(PETSC_COMM_SELF,"out.bin",FILE_MODE_WRITE,&viewer);
1038:   MatView(A,viewer);
1039:   PetscViewerDestroy(&viewer);
1040:   MatDestroy(&A);
1041:   SNESDestroy(&snes_alg);
1042:   VecDestroy(&F_alg);
1043:   MatDestroy(&J);
1044:   MatDestroy(&user.Ybus);
1045:   MatDestroy(&user.Sol);
1046:   VecDestroy(&X);
1047:   VecDestroy(&user.V0);
1048:   DMDestroy(&user.dmgen);
1049:   DMDestroy(&user.dmnet);
1050:   DMDestroy(&user.dmpgrid);
1051:   ISDestroy(&user.is_diff);
1052:   ISDestroy(&user.is_alg);
1053:   TSDestroy(&ts);
1054:   if(user.setisdiff) {
1055:     VecDestroy(&vatol);
1056:   }
1057:   PetscFinalize();
1058:   return(0);
1059: }