Actual source code: ex9bus.c

petsc-3.7.5 2017-01-01
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: {
155:   Userctx        *user;
156:   Vec            X;
157:   PetscScalar    *x,*mat;
158:   PetscInt       idx;
159:   PetscReal      t;

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

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

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

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

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

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

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

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

212:     theta = PETSC_PI/2.0 - delta;

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

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

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

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

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

231:     idx = idx + 4;

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

236:     idx = idx + 2;

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

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

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

250:     idx = idx + 3;
251:   }

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

351:     idx = idx + 9;
352:   }

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

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

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

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

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

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

397:   user->t = t;

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

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

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

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

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

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

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

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

476:     idx = idx + 9;
477:   }

479:   start = user->neqs_gen;

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

488:   MatSeqAIJSetPreallocation(J,0,d_nnz);

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

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


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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

632:     /* Exciter differential equations */

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

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

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

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

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

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

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

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

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

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


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

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

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

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

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

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

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

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

763:   user->t = t;

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

903:   /* Save initial solution */

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

909:   mat[idx] = 0.0;

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

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

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

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

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

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

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

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

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

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

986:   user.alg_flg = PETSC_FALSE;

988:   TSSolve(ts,X);

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

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

1001:   MatZeroEntries(J);

1003:   user.alg_flg = PETSC_TRUE;

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

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

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

1023:   user.alg_flg = PETSC_TRUE;

1025:   TSSolve(ts,X);

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

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