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 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: */
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;
144: PetscScalar Efd,RF,VR,Vr,Vi,Vm;
148: DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
149: DMCompositeScatter(user->dmpgrid,X,Xgen,Xnet);
151: VecGetArrayRead(Xgen,&xgen);
152: VecGetArrayRead(Xnet,&xnet);
154: /* Event for fault-on time */
155: fvalue[0] = t - user->tfaulton;
156: /* Event for fault-off time */
157: fvalue[1] = t - user->tfaultoff;
159: for (i=0; i < ngen; i++) {
160: Efd = xgen[idx+6];
161: RF = xgen[idx+7];
162: VR = xgen[idx+8];
164: Vr = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */
165: Vi = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */
166: Vm = PetscSqrtScalar(Vr*Vr + Vi*Vi);
168: if (!VRatmax[i]) {
169: fvalue[2+2*i] = VRMAX[i] - VR;
170: } else {
171: fvalue[2+2*i] = (VR - KA[i]*RF + KA[i]*KF[i]*Efd/TF[i] - KA[i]*(Vref[i] - Vm))/TA[i];
172: }
173: if (!VRatmin[i]) {
174: fvalue[2+2*i+1] = VRMIN[i] - VR;
175: } else {
176: fvalue[2+2*i+1] = (VR - KA[i]*RF + KA[i]*KF[i]*Efd/TF[i] - KA[i]*(Vref[i] - Vm))/TA[i];
177: }
178: idx = idx+9;
179: }
180: VecRestoreArrayRead(Xgen,&xgen);
181: VecRestoreArrayRead(Xnet,&xnet);
183: DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);
185: return(0);
186: }
188: PetscErrorCode PostEventFunction(TS ts,PetscInt nevents,PetscInt event_list[],PetscReal t,Vec X,PetscBool forwardsolve,void* ctx)
189: {
190: Userctx *user=(Userctx*)ctx;
191: Vec Xgen,Xnet;
192: PetscScalar *xgen,*xnet;
193: PetscInt row_loc,col_loc;
194: PetscScalar val;
196: PetscInt i,idx=0,event_num;
197: PetscScalar fvalue;
198: PetscScalar Efd, RF, VR;
199: PetscScalar Vr,Vi,Vm;
203: DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
204: DMCompositeScatter(user->dmpgrid,X,Xgen,Xnet);
206: VecGetArray(Xgen,&xgen);
207: VecGetArray(Xnet,&xnet);
209: for (i=0; i < nevents; i++) {
210: if (event_list[i] == 0) {
211: /* Apply disturbance - resistive fault at user->faultbus */
212: /* This is done by adding shunt conductance to the diagonal location
213: in the Ybus matrix */
214: row_loc = 2*user->faultbus; col_loc = 2*user->faultbus+1; /* Location for G */
215: val = 1/user->Rfault;
216: MatSetValues(user->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);
217: row_loc = 2*user->faultbus+1; col_loc = 2*user->faultbus; /* Location for G */
218: val = 1/user->Rfault;
219: MatSetValues(user->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);
221: MatAssemblyBegin(user->Ybus,MAT_FINAL_ASSEMBLY);
222: MatAssemblyEnd(user->Ybus,MAT_FINAL_ASSEMBLY);
224: /* Solve the algebraic equations */
225: SNESSolve(user->snes_alg,NULL,X);
226: } else if (event_list[i] == 1) {
227: /* Remove the fault */
228: row_loc = 2*user->faultbus; col_loc = 2*user->faultbus+1;
229: val = -1/user->Rfault;
230: MatSetValues(user->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);
231: row_loc = 2*user->faultbus+1; col_loc = 2*user->faultbus;
232: val = -1/user->Rfault;
233: MatSetValues(user->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);
235: MatAssemblyBegin(user->Ybus,MAT_FINAL_ASSEMBLY);
236: MatAssemblyEnd(user->Ybus,MAT_FINAL_ASSEMBLY);
238: /* Solve the algebraic equations */
239: SNESSolve(user->snes_alg,NULL,X);
241: /* Check the VR derivatives and reset flags if needed */
242: for (i=0; i < ngen; i++) {
243: Efd = xgen[idx+6];
244: RF = xgen[idx+7];
245: VR = xgen[idx+8];
247: Vr = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */
248: Vi = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */
249: Vm = PetscSqrtScalar(Vr*Vr + Vi*Vi);
251: if (VRatmax[i]) {
252: fvalue = (VR - KA[i]*RF + KA[i]*KF[i]*Efd/TF[i] - KA[i]*(Vref[i] - Vm))/TA[i];
253: if (fvalue < 0) {
254: VRatmax[i] = 0;
255: PetscPrintf(PETSC_COMM_SELF,"VR[%d]: dVR_dt went negative on fault clearing at time %g\n",i,t);
256: }
257: }
258: if (VRatmin[i]) {
259: fvalue = (VR - KA[i]*RF + KA[i]*KF[i]*Efd/TF[i] - KA[i]*(Vref[i] - Vm))/TA[i];
261: if (fvalue > 0) {
262: VRatmin[i] = 0;
263: PetscPrintf(PETSC_COMM_SELF,"VR[%d]: dVR_dt went positive on fault clearing at time %g\n",i,t);
264: }
265: }
266: idx = idx+9;
267: }
268: } else {
269: idx = (event_list[i]-2)/2;
270: event_num = (event_list[i]-2)%2;
271: if (event_num == 0) { /* Max VR */
272: if (!VRatmax[idx]) {
273: VRatmax[idx] = 1;
274: PetscPrintf(PETSC_COMM_SELF,"VR[%d]: hit upper limit at time %g\n",idx,t);
275: }
276: else {
277: VRatmax[idx] = 0;
278: PetscPrintf(PETSC_COMM_SELF,"VR[%d]: freeing variable as dVR_dt is negative at time %g\n",idx,t);
279: }
280: } else {
281: if (!VRatmin[idx]) {
282: VRatmin[idx] = 1;
283: PetscPrintf(PETSC_COMM_SELF,"VR[%d]: hit lower limit at time %g\n",idx,t);
284: }
285: else {
286: VRatmin[idx] = 0;
287: PetscPrintf(PETSC_COMM_SELF,"VR[%d]: freeing variable as dVR_dt is positive at time %g\n",idx,t);
288: }
289: }
290: }
291: }
293: VecRestoreArray(Xgen,&xgen);
294: VecRestoreArray(Xnet,&xnet);
296: DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);
298: return(0);
299: }
301: /* Converts from machine frame (dq) to network (phase a real,imag) reference frame */
302: PetscErrorCode dq2ri(PetscScalar Fd,PetscScalar Fq,PetscScalar delta,PetscScalar *Fr, PetscScalar *Fi)
303: {
305: *Fr = Fd*PetscSinScalar(delta) + Fq*PetscCosScalar(delta);
306: *Fi = -Fd*PetscCosScalar(delta) + Fq*PetscSinScalar(delta);
307: return(0);
308: }
310: /* Converts from network frame ([phase a real,imag) to machine (dq) reference frame */
311: PetscErrorCode ri2dq(PetscScalar Fr,PetscScalar Fi,PetscScalar delta,PetscScalar *Fd, PetscScalar *Fq)
312: {
314: *Fd = Fr*PetscSinScalar(delta) - Fi*PetscCosScalar(delta);
315: *Fq = Fr*PetscCosScalar(delta) + Fi*PetscSinScalar(delta);
316: return(0);
317: }
319: /* Saves the solution at each time to a matrix */
320: PetscErrorCode SaveSolution(TS ts)
321: {
322: PetscErrorCode ierr;
323: Userctx *user;
324: Vec X;
325: const PetscScalar *x;
326: PetscScalar *mat;
327: PetscInt idx;
328: PetscReal t;
331: TSGetApplicationContext(ts,&user);
332: TSGetTime(ts,&t);
333: TSGetSolution(ts,&X);
334: idx = user->stepnum*(user->neqs_pgrid+1);
335: MatDenseGetArray(user->Sol,&mat);
336: VecGetArrayRead(X,&x);
337: mat[idx] = t;
338: PetscArraycpy(mat+idx+1,x,user->neqs_pgrid);
339: MatDenseRestoreArray(user->Sol,&mat);
340: VecRestoreArrayRead(X,&x);
341: user->stepnum++;
342: return(0);
343: }
345: PetscErrorCode SetInitialGuess(Vec X,Userctx *user)
346: {
347: PetscErrorCode ierr;
348: Vec Xgen,Xnet;
349: PetscScalar *xgen;
350: const PetscScalar *xnet;
351: PetscInt i,idx=0;
352: PetscScalar Vr,Vi,IGr,IGi,Vm,Vm2;
353: PetscScalar Eqp,Edp,delta;
354: PetscScalar Efd,RF,VR; /* Exciter variables */
355: PetscScalar Id,Iq; /* Generator dq axis currents */
356: PetscScalar theta,Vd,Vq,SE;
359: M[0] = 2*H[0]/w_s; M[1] = 2*H[1]/w_s; M[2] = 2*H[2]/w_s;
360: D[0] = 0.1*M[0]; D[1] = 0.1*M[1]; D[2] = 0.1*M[2];
362: DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
364: /* Network subsystem initialization */
365: VecCopy(user->V0,Xnet);
367: /* Generator subsystem initialization */
368: VecGetArrayWrite(Xgen,&xgen);
369: VecGetArrayRead(Xnet,&xnet);
371: for (i=0; i < ngen; i++) {
372: Vr = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */
373: Vi = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */
374: Vm = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm;
375: IGr = (Vr*PG[i] + Vi*QG[i])/Vm2;
376: IGi = (Vi*PG[i] - Vr*QG[i])/Vm2;
378: delta = PetscAtan2Real(Vi+Xq[i]*IGr,Vr-Xq[i]*IGi); /* Machine angle */
380: theta = PETSC_PI/2.0 - delta;
382: Id = IGr*PetscCosScalar(theta) - IGi*PetscSinScalar(theta); /* d-axis stator current */
383: Iq = IGr*PetscSinScalar(theta) + IGi*PetscCosScalar(theta); /* q-axis stator current */
385: Vd = Vr*PetscCosScalar(theta) - Vi*PetscSinScalar(theta);
386: Vq = Vr*PetscSinScalar(theta) + Vi*PetscCosScalar(theta);
388: Edp = Vd + Rs[i]*Id - Xqp[i]*Iq; /* d-axis transient EMF */
389: Eqp = Vq + Rs[i]*Iq + Xdp[i]*Id; /* q-axis transient EMF */
391: TM[i] = PG[i];
393: /* The generator variables are ordered as [Eqp,Edp,delta,w,Id,Iq] */
394: xgen[idx] = Eqp;
395: xgen[idx+1] = Edp;
396: xgen[idx+2] = delta;
397: xgen[idx+3] = w_s;
399: idx = idx + 4;
401: xgen[idx] = Id;
402: xgen[idx+1] = Iq;
404: idx = idx + 2;
406: /* Exciter */
407: Efd = Eqp + (Xd[i] - Xdp[i])*Id;
408: SE = k1[i]*PetscExpScalar(k2[i]*Efd);
409: VR = KE[i]*Efd + SE;
410: RF = KF[i]*Efd/TF[i];
412: xgen[idx] = Efd;
413: xgen[idx+1] = RF;
414: xgen[idx+2] = VR;
416: Vref[i] = Vm + (VR/KA[i]);
418: VRatmin[i] = VRatmax[i] = 0;
420: idx = idx + 3;
421: }
422: VecRestoreArrayWrite(Xgen,&xgen);
423: VecRestoreArrayRead(Xnet,&xnet);
425: /* VecView(Xgen,0); */
426: DMCompositeGather(user->dmpgrid,INSERT_VALUES,X,Xgen,Xnet);
427: DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);
428: return(0);
429: }
431: /* Computes F = [f(x,y);g(x,y)] */
432: PetscErrorCode ResidualFunction(Vec X, Vec F, Userctx *user)
433: {
434: PetscErrorCode ierr;
435: Vec Xgen,Xnet,Fgen,Fnet;
436: const PetscScalar *xgen,*xnet;
437: PetscScalar *fgen,*fnet;
438: PetscInt i,idx=0;
439: PetscScalar Vr,Vi,Vm,Vm2;
440: PetscScalar Eqp,Edp,delta,w; /* Generator variables */
441: PetscScalar Efd,RF,VR; /* Exciter variables */
442: PetscScalar Id,Iq; /* Generator dq axis currents */
443: PetscScalar Vd,Vq,SE;
444: PetscScalar IGr,IGi,IDr,IDi;
445: PetscScalar Zdq_inv[4],det;
446: PetscScalar PD,QD,Vm0,*v0;
447: PetscInt k;
450: VecZeroEntries(F);
451: DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
452: DMCompositeGetLocalVectors(user->dmpgrid,&Fgen,&Fnet);
453: DMCompositeScatter(user->dmpgrid,X,Xgen,Xnet);
454: DMCompositeScatter(user->dmpgrid,F,Fgen,Fnet);
456: /* Network current balance residual IG + Y*V + IL = 0. Only YV is added here.
457: The generator current injection, IG, and load current injection, ID are added later
458: */
459: /* Note that the values in Ybus are stored assuming the imaginary current balance
460: equation is ordered first followed by real current balance equation for each bus.
461: Thus imaginary current contribution goes in location 2*i, and
462: real current contribution in 2*i+1
463: */
464: MatMult(user->Ybus,Xnet,Fnet);
466: VecGetArrayRead(Xgen,&xgen);
467: VecGetArrayRead(Xnet,&xnet);
468: VecGetArrayWrite(Fgen,&fgen);
469: VecGetArrayWrite(Fnet,&fnet);
471: /* Generator subsystem */
472: for (i=0; i < ngen; i++) {
473: Eqp = xgen[idx];
474: Edp = xgen[idx+1];
475: delta = xgen[idx+2];
476: w = xgen[idx+3];
477: Id = xgen[idx+4];
478: Iq = xgen[idx+5];
479: Efd = xgen[idx+6];
480: RF = xgen[idx+7];
481: VR = xgen[idx+8];
483: /* Generator differential equations */
484: fgen[idx] = (-Eqp - (Xd[i] - Xdp[i])*Id + Efd)/Td0p[i];
485: fgen[idx+1] = (-Edp + (Xq[i] - Xqp[i])*Iq)/Tq0p[i];
486: fgen[idx+2] = w - w_s;
487: fgen[idx+3] = (TM[i] - Edp*Id - Eqp*Iq - (Xqp[i] - Xdp[i])*Id*Iq - D[i]*(w - w_s))/M[i];
489: Vr = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */
490: Vi = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */
492: ri2dq(Vr,Vi,delta,&Vd,&Vq);
493: /* Algebraic equations for stator currents */
494: det = Rs[i]*Rs[i] + Xdp[i]*Xqp[i];
496: Zdq_inv[0] = Rs[i]/det;
497: Zdq_inv[1] = Xqp[i]/det;
498: Zdq_inv[2] = -Xdp[i]/det;
499: Zdq_inv[3] = Rs[i]/det;
501: fgen[idx+4] = Zdq_inv[0]*(-Edp + Vd) + Zdq_inv[1]*(-Eqp + Vq) + Id;
502: fgen[idx+5] = Zdq_inv[2]*(-Edp + Vd) + Zdq_inv[3]*(-Eqp + Vq) + Iq;
504: /* Add generator current injection to network */
505: dq2ri(Id,Iq,delta,&IGr,&IGi);
507: fnet[2*gbus[i]] -= IGi;
508: fnet[2*gbus[i]+1] -= IGr;
510: Vm = PetscSqrtScalar(Vd*Vd + Vq*Vq);
512: SE = k1[i]*PetscExpScalar(k2[i]*Efd);
514: /* Exciter differential equations */
515: fgen[idx+6] = (-KE[i]*Efd - SE + VR)/TE[i];
516: fgen[idx+7] = (-RF + KF[i]*Efd/TF[i])/TF[i];
517: if (VRatmax[i]) fgen[idx+8] = VR - VRMAX[i];
518: else if (VRatmin[i]) fgen[idx+8] = VRMIN[i] - VR;
519: else fgen[idx+8] = (-VR + KA[i]*RF - KA[i]*KF[i]*Efd/TF[i] + KA[i]*(Vref[i] - Vm))/TA[i];
521: idx = idx + 9;
522: }
524: VecGetArray(user->V0,&v0);
525: for (i=0; i < nload; i++) {
526: Vr = xnet[2*lbus[i]]; /* Real part of load bus voltage */
527: Vi = xnet[2*lbus[i]+1]; /* Imaginary part of the load bus voltage */
528: Vm = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm;
529: Vm0 = PetscSqrtScalar(v0[2*lbus[i]]*v0[2*lbus[i]] + v0[2*lbus[i]+1]*v0[2*lbus[i]+1]);
530: PD = QD = 0.0;
531: for (k=0; k < ld_nsegsp[i]; k++) PD += ld_alphap[k]*PD0[i]*PetscPowScalar((Vm/Vm0),ld_betap[k]);
532: for (k=0; k < ld_nsegsq[i]; k++) QD += ld_alphaq[k]*QD0[i]*PetscPowScalar((Vm/Vm0),ld_betaq[k]);
534: /* Load currents */
535: IDr = (PD*Vr + QD*Vi)/Vm2;
536: IDi = (-QD*Vr + PD*Vi)/Vm2;
538: fnet[2*lbus[i]] += IDi;
539: fnet[2*lbus[i]+1] += IDr;
540: }
541: VecRestoreArray(user->V0,&v0);
543: VecRestoreArrayRead(Xgen,&xgen);
544: VecRestoreArrayRead(Xnet,&xnet);
545: VecRestoreArrayWrite(Fgen,&fgen);
546: VecRestoreArrayWrite(Fnet,&fnet);
548: DMCompositeGather(user->dmpgrid,INSERT_VALUES,F,Fgen,Fnet);
549: DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);
550: DMCompositeRestoreLocalVectors(user->dmpgrid,&Fgen,&Fnet);
551: return(0);
552: }
554: /* f(x,y)
555: g(x,y)
556: */
557: PetscErrorCode RHSFunction(TS ts,PetscReal t, Vec X, Vec F, void *ctx)
558: {
560: Userctx *user=(Userctx*)ctx;
563: user->t = t;
564: ResidualFunction(X,F,user);
565: return(0);
566: }
568: /* f(x,y) - \dot{x}
569: g(x,y) = 0
570: */
571: PetscErrorCode IFunction(TS ts,PetscReal t, Vec X, Vec Xdot, Vec F, void *ctx)
572: {
573: PetscErrorCode ierr;
574: PetscScalar *f;
575: const PetscScalar *xdot;
576: PetscInt i;
579: RHSFunction(ts,t,X,F,ctx);
580: VecScale(F,-1.0);
581: VecGetArray(F,&f);
582: VecGetArrayRead(Xdot,&xdot);
583: for (i=0;i < ngen;i++) {
584: f[9*i] += xdot[9*i];
585: f[9*i+1] += xdot[9*i+1];
586: f[9*i+2] += xdot[9*i+2];
587: f[9*i+3] += xdot[9*i+3];
588: f[9*i+6] += xdot[9*i+6];
589: f[9*i+7] += xdot[9*i+7];
590: f[9*i+8] += xdot[9*i+8];
591: }
592: VecRestoreArray(F,&f);
593: VecRestoreArrayRead(Xdot,&xdot);
594: return(0);
595: }
597: /* This function is used for solving the algebraic system only during fault on and
598: off times. It computes the entire F and then zeros out the part corresponding to
599: differential equations
600: F = [0;g(y)];
601: */
602: PetscErrorCode AlgFunction(SNES snes, Vec X, Vec F, void *ctx)
603: {
605: Userctx *user=(Userctx*)ctx;
606: PetscScalar *f;
607: PetscInt i;
610: ResidualFunction(X,F,user);
611: VecGetArray(F,&f);
612: for (i=0; i < ngen; i++) {
613: f[9*i] = 0;
614: f[9*i+1] = 0;
615: f[9*i+2] = 0;
616: f[9*i+3] = 0;
617: f[9*i+6] = 0;
618: f[9*i+7] = 0;
619: f[9*i+8] = 0;
620: }
621: VecRestoreArray(F,&f);
622: return(0);
623: }
625: PetscErrorCode PostStage(TS ts, PetscReal t, PetscInt i, Vec *X)
626: {
628: Userctx *user;
631: TSGetApplicationContext(ts,&user);
632: SNESSolve(user->snes_alg,NULL,X[i]);
633: return(0);
634: }
636: PetscErrorCode PostEvaluate(TS ts)
637: {
639: Userctx *user;
640: Vec X;
643: TSGetApplicationContext(ts,&user);
644: TSGetSolution(ts,&X);
645: SNESSolve(user->snes_alg,NULL,X);
646: return(0);
647: }
650: PetscErrorCode PreallocateJacobian(Mat J, Userctx *user)
651: {
653: PetscInt *d_nnz;
654: PetscInt i,idx=0,start=0;
655: PetscInt ncols;
658: PetscMalloc1(user->neqs_pgrid,&d_nnz);
659: for (i=0; i<user->neqs_pgrid; i++) d_nnz[i] = 0;
660: /* Generator subsystem */
661: for (i=0; i < ngen; i++) {
663: d_nnz[idx] += 3;
664: d_nnz[idx+1] += 2;
665: d_nnz[idx+2] += 2;
666: d_nnz[idx+3] += 5;
667: d_nnz[idx+4] += 6;
668: d_nnz[idx+5] += 6;
670: d_nnz[user->neqs_gen+2*gbus[i]] += 3;
671: d_nnz[user->neqs_gen+2*gbus[i]+1] += 3;
673: d_nnz[idx+6] += 2;
674: d_nnz[idx+7] += 2;
675: d_nnz[idx+8] += 5;
677: idx = idx + 9;
678: }
679: start = user->neqs_gen;
681: for (i=0; i < nbus; i++) {
682: MatGetRow(user->Ybus,2*i,&ncols,NULL,NULL);
683: d_nnz[start+2*i] += ncols;
684: d_nnz[start+2*i+1] += ncols;
685: MatRestoreRow(user->Ybus,2*i,&ncols,NULL,NULL);
686: }
687: MatSeqAIJSetPreallocation(J,0,d_nnz);
688: PetscFree(d_nnz);
689: return(0);
690: }
692: /*
693: J = [df_dx, df_dy
694: dg_dx, dg_dy]
695: */
696: PetscErrorCode ResidualJacobian(Vec X,Mat J,Mat B,void *ctx)
697: {
698: PetscErrorCode ierr;
699: Userctx *user = (Userctx*)ctx;
700: Vec Xgen,Xnet;
701: const PetscScalar *xgen,*xnet;
702: PetscInt i,idx=0;
703: PetscScalar Vr,Vi,Vm,Vm2;
704: PetscScalar Eqp,Edp,delta; /* Generator variables */
705: PetscScalar Efd;
706: PetscScalar Id,Iq; /* Generator dq axis currents */
707: PetscScalar Vd,Vq;
708: PetscScalar val[10];
709: PetscInt row[2],col[10];
710: PetscInt net_start = user->neqs_gen;
711: PetscScalar Zdq_inv[4],det;
712: PetscScalar dVd_dVr,dVd_dVi,dVq_dVr,dVq_dVi,dVd_ddelta,dVq_ddelta;
713: PetscScalar dIGr_ddelta,dIGi_ddelta,dIGr_dId,dIGr_dIq,dIGi_dId,dIGi_dIq;
714: PetscScalar dSE_dEfd;
715: PetscScalar dVm_dVd,dVm_dVq,dVm_dVr,dVm_dVi;
716: PetscInt ncols;
717: const PetscInt *cols;
718: const PetscScalar *yvals;
719: PetscInt k;
720: PetscScalar PD,QD,Vm0,Vm4;
721: const PetscScalar *v0;
722: PetscScalar dPD_dVr,dPD_dVi,dQD_dVr,dQD_dVi;
723: PetscScalar dIDr_dVr,dIDr_dVi,dIDi_dVr,dIDi_dVi;
727: MatZeroEntries(B);
728: DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
729: DMCompositeScatter(user->dmpgrid,X,Xgen,Xnet);
731: VecGetArrayRead(Xgen,&xgen);
732: VecGetArrayRead(Xnet,&xnet);
734: /* Generator subsystem */
735: for (i=0; i < ngen; i++) {
736: Eqp = xgen[idx];
737: Edp = xgen[idx+1];
738: delta = xgen[idx+2];
739: Id = xgen[idx+4];
740: Iq = xgen[idx+5];
741: Efd = xgen[idx+6];
743: /* fgen[idx] = (-Eqp - (Xd[i] - Xdp[i])*Id + Efd)/Td0p[i]; */
744: row[0] = idx;
745: col[0] = idx; col[1] = idx+4; col[2] = idx+6;
746: val[0] = -1/ Td0p[i]; val[1] = -(Xd[i] - Xdp[i])/ Td0p[i]; val[2] = 1/Td0p[i];
748: MatSetValues(J,1,row,3,col,val,INSERT_VALUES);
750: /* fgen[idx+1] = (-Edp + (Xq[i] - Xqp[i])*Iq)/Tq0p[i]; */
751: row[0] = idx + 1;
752: col[0] = idx + 1; col[1] = idx+5;
753: val[0] = -1/Tq0p[i]; val[1] = (Xq[i] - Xqp[i])/Tq0p[i];
754: MatSetValues(J,1,row,2,col,val,INSERT_VALUES);
756: /* fgen[idx+2] = w - w_s; */
757: row[0] = idx + 2;
758: col[0] = idx + 2; col[1] = idx + 3;
759: val[0] = 0; val[1] = 1;
760: MatSetValues(J,1,row,2,col,val,INSERT_VALUES);
762: /* fgen[idx+3] = (TM[i] - Edp*Id - Eqp*Iq - (Xqp[i] - Xdp[i])*Id*Iq - D[i]*(w - w_s))/M[i]; */
763: row[0] = idx + 3;
764: col[0] = idx; col[1] = idx + 1; col[2] = idx + 3; col[3] = idx + 4; col[4] = idx + 5;
765: 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];
766: MatSetValues(J,1,row,5,col,val,INSERT_VALUES);
768: Vr = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */
769: Vi = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */
770: ri2dq(Vr,Vi,delta,&Vd,&Vq);
772: det = Rs[i]*Rs[i] + Xdp[i]*Xqp[i];
774: Zdq_inv[0] = Rs[i]/det;
775: Zdq_inv[1] = Xqp[i]/det;
776: Zdq_inv[2] = -Xdp[i]/det;
777: Zdq_inv[3] = Rs[i]/det;
779: dVd_dVr = PetscSinScalar(delta); dVd_dVi = -PetscCosScalar(delta);
780: dVq_dVr = PetscCosScalar(delta); dVq_dVi = PetscSinScalar(delta);
781: dVd_ddelta = Vr*PetscCosScalar(delta) + Vi*PetscSinScalar(delta);
782: dVq_ddelta = -Vr*PetscSinScalar(delta) + Vi*PetscCosScalar(delta);
784: /* fgen[idx+4] = Zdq_inv[0]*(-Edp + Vd) + Zdq_inv[1]*(-Eqp + Vq) + Id; */
785: row[0] = idx+4;
786: col[0] = idx; col[1] = idx+1; col[2] = idx + 2;
787: val[0] = -Zdq_inv[1]; val[1] = -Zdq_inv[0]; val[2] = Zdq_inv[0]*dVd_ddelta + Zdq_inv[1]*dVq_ddelta;
788: col[3] = idx + 4; col[4] = net_start+2*gbus[i]; col[5] = net_start + 2*gbus[i]+1;
789: 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;
790: MatSetValues(J,1,row,6,col,val,INSERT_VALUES);
792: /* fgen[idx+5] = Zdq_inv[2]*(-Edp + Vd) + Zdq_inv[3]*(-Eqp + Vq) + Iq; */
793: row[0] = idx+5;
794: col[0] = idx; col[1] = idx+1; col[2] = idx + 2;
795: val[0] = -Zdq_inv[3]; val[1] = -Zdq_inv[2]; val[2] = Zdq_inv[2]*dVd_ddelta + Zdq_inv[3]*dVq_ddelta;
796: col[3] = idx + 5; col[4] = net_start+2*gbus[i]; col[5] = net_start + 2*gbus[i]+1;
797: 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;
798: MatSetValues(J,1,row,6,col,val,INSERT_VALUES);
800: dIGr_ddelta = Id*PetscCosScalar(delta) - Iq*PetscSinScalar(delta);
801: dIGi_ddelta = Id*PetscSinScalar(delta) + Iq*PetscCosScalar(delta);
802: dIGr_dId = PetscSinScalar(delta); dIGr_dIq = PetscCosScalar(delta);
803: dIGi_dId = -PetscCosScalar(delta); dIGi_dIq = PetscSinScalar(delta);
805: /* fnet[2*gbus[i]] -= IGi; */
806: row[0] = net_start + 2*gbus[i];
807: col[0] = idx+2; col[1] = idx + 4; col[2] = idx + 5;
808: val[0] = -dIGi_ddelta; val[1] = -dIGi_dId; val[2] = -dIGi_dIq;
809: MatSetValues(J,1,row,3,col,val,INSERT_VALUES);
811: /* fnet[2*gbus[i]+1] -= IGr; */
812: row[0] = net_start + 2*gbus[i]+1;
813: col[0] = idx+2; col[1] = idx + 4; col[2] = idx + 5;
814: val[0] = -dIGr_ddelta; val[1] = -dIGr_dId; val[2] = -dIGr_dIq;
815: MatSetValues(J,1,row,3,col,val,INSERT_VALUES);
817: Vm = PetscSqrtScalar(Vd*Vd + Vq*Vq);
819: /* fgen[idx+6] = (-KE[i]*Efd - SE + VR)/TE[i]; */
820: /* SE = k1[i]*PetscExpScalar(k2[i]*Efd); */
822: dSE_dEfd = k1[i]*k2[i]*PetscExpScalar(k2[i]*Efd);
824: row[0] = idx + 6;
825: col[0] = idx + 6; col[1] = idx + 8;
826: val[0] = (-KE[i] - dSE_dEfd)/TE[i]; val[1] = 1/TE[i];
827: MatSetValues(J,1,row,2,col,val,INSERT_VALUES);
829: /* Exciter differential equations */
831: /* fgen[idx+7] = (-RF + KF[i]*Efd/TF[i])/TF[i]; */
832: row[0] = idx + 7;
833: col[0] = idx + 6; col[1] = idx + 7;
834: val[0] = (KF[i]/TF[i])/TF[i]; val[1] = -1/TF[i];
835: MatSetValues(J,1,row,2,col,val,INSERT_VALUES);
837: /* fgen[idx+8] = (-VR + KA[i]*RF - KA[i]*KF[i]*Efd/TF[i] + KA[i]*(Vref[i] - Vm))/TA[i]; */
838: /* Vm = (Vd^2 + Vq^2)^0.5; */
840: row[0] = idx + 8;
841: if (VRatmax[i]) {
842: col[0] = idx + 8; val[0] = 1.0;
843: MatSetValues(J,1,row,1,col,val,INSERT_VALUES);
844: } else if (VRatmin[i]) {
845: col[0] = idx + 8; val[0] = -1.0;
846: MatSetValues(J,1,row,1,col,val,INSERT_VALUES);
847: } else {
848: dVm_dVd = Vd/Vm; dVm_dVq = Vq/Vm;
849: dVm_dVr = dVm_dVd*dVd_dVr + dVm_dVq*dVq_dVr;
850: dVm_dVi = dVm_dVd*dVd_dVi + dVm_dVq*dVq_dVi;
851: row[0] = idx + 8;
852: col[0] = idx + 6; col[1] = idx + 7; col[2] = idx + 8;
853: val[0] = -(KA[i]*KF[i]/TF[i])/TA[i]; val[1] = KA[i]/TA[i]; val[2] = -1/TA[i];
854: col[3] = net_start + 2*gbus[i]; col[4] = net_start + 2*gbus[i]+1;
855: val[3] = -KA[i]*dVm_dVr/TA[i]; val[4] = -KA[i]*dVm_dVi/TA[i];
856: MatSetValues(J,1,row,5,col,val,INSERT_VALUES);
857: }
858: idx = idx + 9;
859: }
861: for (i=0; i<nbus; i++) {
862: MatGetRow(user->Ybus,2*i,&ncols,&cols,&yvals);
863: row[0] = net_start + 2*i;
864: for (k=0; k<ncols; k++) {
865: col[k] = net_start + cols[k];
866: val[k] = yvals[k];
867: }
868: MatSetValues(J,1,row,ncols,col,val,INSERT_VALUES);
869: MatRestoreRow(user->Ybus,2*i,&ncols,&cols,&yvals);
871: MatGetRow(user->Ybus,2*i+1,&ncols,&cols,&yvals);
872: row[0] = net_start + 2*i+1;
873: for (k=0; k<ncols; k++) {
874: col[k] = net_start + cols[k];
875: val[k] = yvals[k];
876: }
877: MatSetValues(J,1,row,ncols,col,val,INSERT_VALUES);
878: MatRestoreRow(user->Ybus,2*i+1,&ncols,&cols,&yvals);
879: }
881: MatAssemblyBegin(J,MAT_FLUSH_ASSEMBLY);
882: MatAssemblyEnd(J,MAT_FLUSH_ASSEMBLY);
884: VecGetArrayRead(user->V0,&v0);
885: for (i=0; i < nload; i++) {
886: Vr = xnet[2*lbus[i]]; /* Real part of load bus voltage */
887: Vi = xnet[2*lbus[i]+1]; /* Imaginary part of the load bus voltage */
888: Vm = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm; Vm4 = Vm2*Vm2;
889: Vm0 = PetscSqrtScalar(v0[2*lbus[i]]*v0[2*lbus[i]] + v0[2*lbus[i]+1]*v0[2*lbus[i]+1]);
890: PD = QD = 0.0;
891: dPD_dVr = dPD_dVi = dQD_dVr = dQD_dVi = 0.0;
892: for (k=0; k < ld_nsegsp[i]; k++) {
893: PD += ld_alphap[k]*PD0[i]*PetscPowScalar((Vm/Vm0),ld_betap[k]);
894: dPD_dVr += ld_alphap[k]*ld_betap[k]*PD0[i]*PetscPowScalar((1/Vm0),ld_betap[k])*Vr*PetscPowScalar(Vm,(ld_betap[k]-2));
895: dPD_dVi += ld_alphap[k]*ld_betap[k]*PD0[i]*PetscPowScalar((1/Vm0),ld_betap[k])*Vi*PetscPowScalar(Vm,(ld_betap[k]-2));
896: }
897: for (k=0; k < ld_nsegsq[i]; k++) {
898: QD += ld_alphaq[k]*QD0[i]*PetscPowScalar((Vm/Vm0),ld_betaq[k]);
899: dQD_dVr += ld_alphaq[k]*ld_betaq[k]*QD0[i]*PetscPowScalar((1/Vm0),ld_betaq[k])*Vr*PetscPowScalar(Vm,(ld_betaq[k]-2));
900: dQD_dVi += ld_alphaq[k]*ld_betaq[k]*QD0[i]*PetscPowScalar((1/Vm0),ld_betaq[k])*Vi*PetscPowScalar(Vm,(ld_betaq[k]-2));
901: }
903: /* IDr = (PD*Vr + QD*Vi)/Vm2; */
904: /* IDi = (-QD*Vr + PD*Vi)/Vm2; */
906: dIDr_dVr = (dPD_dVr*Vr + dQD_dVr*Vi + PD)/Vm2 - ((PD*Vr + QD*Vi)*2*Vr)/Vm4;
907: dIDr_dVi = (dPD_dVi*Vr + dQD_dVi*Vi + QD)/Vm2 - ((PD*Vr + QD*Vi)*2*Vi)/Vm4;
909: dIDi_dVr = (-dQD_dVr*Vr + dPD_dVr*Vi - QD)/Vm2 - ((-QD*Vr + PD*Vi)*2*Vr)/Vm4;
910: dIDi_dVi = (-dQD_dVi*Vr + dPD_dVi*Vi + PD)/Vm2 - ((-QD*Vr + PD*Vi)*2*Vi)/Vm4;
913: /* fnet[2*lbus[i]] += IDi; */
914: row[0] = net_start + 2*lbus[i];
915: col[0] = net_start + 2*lbus[i]; col[1] = net_start + 2*lbus[i]+1;
916: val[0] = dIDi_dVr; val[1] = dIDi_dVi;
917: MatSetValues(J,1,row,2,col,val,ADD_VALUES);
918: /* fnet[2*lbus[i]+1] += IDr; */
919: row[0] = net_start + 2*lbus[i]+1;
920: col[0] = net_start + 2*lbus[i]; col[1] = net_start + 2*lbus[i]+1;
921: val[0] = dIDr_dVr; val[1] = dIDr_dVi;
922: MatSetValues(J,1,row,2,col,val,ADD_VALUES);
923: }
924: VecRestoreArrayRead(user->V0,&v0);
926: VecRestoreArrayRead(Xgen,&xgen);
927: VecRestoreArrayRead(Xnet,&xnet);
929: DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);
931: MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
932: MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
933: return(0);
934: }
936: /*
937: J = [I, 0
938: dg_dx, dg_dy]
939: */
940: PetscErrorCode AlgJacobian(SNES snes,Vec X,Mat A,Mat B,void *ctx)
941: {
943: Userctx *user=(Userctx*)ctx;
946: ResidualJacobian(X,A,B,ctx);
947: MatSetOption(A,MAT_KEEP_NONZERO_PATTERN,PETSC_TRUE);
948: MatZeroRowsIS(A,user->is_diff,1.0,NULL,NULL);
949: return(0);
950: }
952: /*
953: J = [-df_dx, -df_dy
954: dg_dx, dg_dy]
955: */
957: PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec X,Mat A,Mat B,void *ctx)
958: {
960: Userctx *user=(Userctx*)ctx;
963: user->t = t;
965: ResidualJacobian(X,A,B,user);
967: return(0);
968: }
970: /*
971: J = [df_dx-aI, df_dy
972: dg_dx, dg_dy]
973: */
975: PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,Userctx *user)
976: {
978: PetscScalar atmp = (PetscScalar) a;
979: PetscInt i,row;
982: user->t = t;
984: RHSJacobian(ts,t,X,A,B,user);
985: MatScale(B,-1.0);
986: for (i=0;i < ngen;i++) {
987: row = 9*i;
988: MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
989: row = 9*i+1;
990: MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
991: row = 9*i+2;
992: MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
993: row = 9*i+3;
994: MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
995: row = 9*i+6;
996: MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
997: row = 9*i+7;
998: MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
999: row = 9*i+8;
1000: MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
1001: }
1002: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
1003: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
1004: return(0);
1005: }
1007: int main(int argc,char **argv)
1008: {
1009: TS ts;
1010: SNES snes_alg;
1011: PetscErrorCode ierr;
1012: PetscMPIInt size;
1013: Userctx user;
1014: PetscViewer Xview,Ybusview,viewer;
1015: Vec X,F_alg;
1016: Mat J,A;
1017: PetscInt i,idx,*idx2;
1018: Vec Xdot;
1019: PetscScalar *x,*mat,*amat;
1020: const PetscScalar *rmat;
1021: Vec vatol;
1022: PetscInt *direction;
1023: PetscBool *terminate;
1024: const PetscInt *idx3;
1025: PetscScalar *vatoli;
1026: PetscInt k;
1029: PetscInitialize(&argc,&argv,"petscoptions",help);if (ierr) return ierr;
1030: MPI_Comm_size(PETSC_COMM_WORLD,&size);
1031: if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs");
1033: user.neqs_gen = 9*ngen; /* # eqs. for generator subsystem */
1034: user.neqs_net = 2*nbus; /* # eqs. for network subsystem */
1035: user.neqs_pgrid = user.neqs_gen + user.neqs_net;
1037: /* Create indices for differential and algebraic equations */
1039: PetscMalloc1(7*ngen,&idx2);
1040: for (i=0; i<ngen; i++) {
1041: 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;
1042: idx2[7*i+4] = 9*i+6; idx2[7*i+5] = 9*i+7; idx2[7*i+6] = 9*i+8;
1043: }
1044: ISCreateGeneral(PETSC_COMM_WORLD,7*ngen,idx2,PETSC_COPY_VALUES,&user.is_diff);
1045: ISComplement(user.is_diff,0,user.neqs_pgrid,&user.is_alg);
1046: PetscFree(idx2);
1048: /* Read initial voltage vector and Ybus */
1049: PetscViewerBinaryOpen(PETSC_COMM_WORLD,"X.bin",FILE_MODE_READ,&Xview);
1050: PetscViewerBinaryOpen(PETSC_COMM_WORLD,"Ybus.bin",FILE_MODE_READ,&Ybusview);
1052: VecCreate(PETSC_COMM_WORLD,&user.V0);
1053: VecSetSizes(user.V0,PETSC_DECIDE,user.neqs_net);
1054: VecLoad(user.V0,Xview);
1056: MatCreate(PETSC_COMM_WORLD,&user.Ybus);
1057: MatSetSizes(user.Ybus,PETSC_DECIDE,PETSC_DECIDE,user.neqs_net,user.neqs_net);
1058: MatSetType(user.Ybus,MATBAIJ);
1059: /* MatSetBlockSize(user.Ybus,2); */
1060: MatLoad(user.Ybus,Ybusview);
1062: /* Set run time options */
1063: PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Transient stability fault options","");
1064: {
1065: user.tfaulton = 1.0;
1066: user.tfaultoff = 1.2;
1067: user.Rfault = 0.0001;
1068: user.setisdiff = PETSC_FALSE;
1069: user.semiexplicit = PETSC_FALSE;
1070: user.faultbus = 8;
1071: PetscOptionsReal("-tfaulton","","",user.tfaulton,&user.tfaulton,NULL);
1072: PetscOptionsReal("-tfaultoff","","",user.tfaultoff,&user.tfaultoff,NULL);
1073: PetscOptionsInt("-faultbus","","",user.faultbus,&user.faultbus,NULL);
1074: user.t0 = 0.0;
1075: user.tmax = 5.0;
1076: PetscOptionsReal("-t0","","",user.t0,&user.t0,NULL);
1077: PetscOptionsReal("-tmax","","",user.tmax,&user.tmax,NULL);
1078: PetscOptionsBool("-setisdiff","","",user.setisdiff,&user.setisdiff,NULL);
1079: PetscOptionsBool("-dae_semiexplicit","","",user.semiexplicit,&user.semiexplicit,NULL);
1080: }
1081: PetscOptionsEnd();
1083: PetscViewerDestroy(&Xview);
1084: PetscViewerDestroy(&Ybusview);
1086: /* Create DMs for generator and network subsystems */
1087: DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,user.neqs_gen,1,1,NULL,&user.dmgen);
1088: DMSetOptionsPrefix(user.dmgen,"dmgen_");
1089: DMSetFromOptions(user.dmgen);
1090: DMSetUp(user.dmgen);
1091: DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,user.neqs_net,1,1,NULL,&user.dmnet);
1092: DMSetOptionsPrefix(user.dmnet,"dmnet_");
1093: DMSetFromOptions(user.dmnet);
1094: DMSetUp(user.dmnet);
1095: /* Create a composite DM packer and add the two DMs */
1096: DMCompositeCreate(PETSC_COMM_WORLD,&user.dmpgrid);
1097: DMSetOptionsPrefix(user.dmpgrid,"pgrid_");
1098: DMCompositeAddDM(user.dmpgrid,user.dmgen);
1099: DMCompositeAddDM(user.dmpgrid,user.dmnet);
1101: DMCreateGlobalVector(user.dmpgrid,&X);
1103: MatCreate(PETSC_COMM_WORLD,&J);
1104: MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,user.neqs_pgrid,user.neqs_pgrid);
1105: MatSetFromOptions(J);
1106: PreallocateJacobian(J,&user);
1108: /* Create matrix to save solutions at each time step */
1109: user.stepnum = 0;
1111: MatCreateSeqDense(PETSC_COMM_SELF,user.neqs_pgrid+1,1002,NULL,&user.Sol);
1112: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1113: Create timestepping solver context
1114: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1115: TSCreate(PETSC_COMM_WORLD,&ts);
1116: TSSetProblemType(ts,TS_NONLINEAR);
1117: if (user.semiexplicit) {
1118: TSSetType(ts,TSRK);
1119: TSSetRHSFunction(ts,NULL,RHSFunction,&user);
1120: TSSetRHSJacobian(ts,J,J,RHSJacobian,&user);
1121: } else {
1122: TSSetType(ts,TSCN);
1123: TSSetEquationType(ts,TS_EQ_DAE_IMPLICIT_INDEX1);
1124: TSSetIFunction(ts,NULL,(TSIFunction) IFunction,&user);
1125: TSSetIJacobian(ts,J,J,(TSIJacobian)IJacobian,&user);
1126: }
1127: TSSetApplicationContext(ts,&user);
1129: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1130: Set initial conditions
1131: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1132: SetInitialGuess(X,&user);
1133: /* Just to set up the Jacobian structure */
1135: VecDuplicate(X,&Xdot);
1136: IJacobian(ts,0.0,X,Xdot,0.0,J,J,&user);
1137: VecDestroy(&Xdot);
1139: /* Save initial solution */
1141: idx=user.stepnum*(user.neqs_pgrid+1);
1142: MatDenseGetArray(user.Sol,&mat);
1143: VecGetArray(X,&x);
1145: mat[idx] = 0.0;
1147: PetscArraycpy(mat+idx+1,x,user.neqs_pgrid);
1148: MatDenseRestoreArray(user.Sol,&mat);
1149: VecRestoreArray(X,&x);
1150: user.stepnum++;
1152: TSSetMaxTime(ts,user.tmax);
1153: TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);
1154: TSSetTimeStep(ts,0.01);
1155: TSSetFromOptions(ts);
1156: TSSetPostStep(ts,SaveSolution);
1157: TSSetSolution(ts,X);
1159: PetscMalloc1((2*ngen+2),&direction);
1160: PetscMalloc1((2*ngen+2),&terminate);
1161: direction[0] = direction[1] = 1;
1162: terminate[0] = terminate[1] = PETSC_FALSE;
1163: for (i=0; i < ngen;i++) {
1164: direction[2+2*i] = -1; direction[2+2*i+1] = 1;
1165: terminate[2+2*i] = terminate[2+2*i+1] = PETSC_FALSE;
1166: }
1168: TSSetEventHandler(ts,2*ngen+2,direction,terminate,EventFunction,PostEventFunction,(void*)&user);
1170: if (user.semiexplicit) {
1171: /* Use a semi-explicit approach with the time-stepping done by an explicit method and the
1172: algrebraic part solved via PostStage and PostEvaluate callbacks
1173: */
1174: TSSetType(ts,TSRK);
1175: TSSetPostStage(ts,PostStage);
1176: TSSetPostEvaluate(ts,PostEvaluate);
1177: }
1180: if (user.setisdiff) {
1181: /* Create vector of absolute tolerances and set the algebraic part to infinity */
1182: VecDuplicate(X,&vatol);
1183: VecSet(vatol,100000.0);
1184: VecGetArray(vatol,&vatoli);
1185: ISGetIndices(user.is_diff,&idx3);
1186: for (k=0; k < 7*ngen; k++) vatoli[idx3[k]] = 1e-2;
1187: VecRestoreArray(vatol,&vatoli);
1188: }
1190: /* Create the nonlinear solver for solving the algebraic system */
1191: /* Note that although the algebraic system needs to be solved only for
1192: Idq and V, we reuse the entire system including xgen. The xgen
1193: variables are held constant by setting their residuals to 0 and
1194: putting a 1 on the Jacobian diagonal for xgen rows
1195: */
1197: VecDuplicate(X,&F_alg);
1198: SNESCreate(PETSC_COMM_WORLD,&snes_alg);
1199: SNESSetFunction(snes_alg,F_alg,AlgFunction,&user);
1200: SNESSetJacobian(snes_alg,J,J,AlgJacobian,&user);
1201: SNESSetFromOptions(snes_alg);
1203: user.snes_alg=snes_alg;
1205: /* Solve */
1206: TSSolve(ts,X);
1208: MatAssemblyBegin(user.Sol,MAT_FINAL_ASSEMBLY);
1209: MatAssemblyEnd(user.Sol,MAT_FINAL_ASSEMBLY);
1211: MatCreateSeqDense(PETSC_COMM_SELF,user.neqs_pgrid+1,user.stepnum,NULL,&A);
1212: MatDenseGetArrayRead(user.Sol,&rmat);
1213: MatDenseGetArray(A,&amat);
1214: PetscArraycpy(amat,rmat,user.stepnum*(user.neqs_pgrid+1));
1215: MatDenseRestoreArray(A,&amat);
1216: MatDenseRestoreArrayRead(user.Sol,&rmat);
1217: PetscViewerBinaryOpen(PETSC_COMM_SELF,"out.bin",FILE_MODE_WRITE,&viewer);
1218: MatView(A,viewer);
1219: PetscViewerDestroy(&viewer);
1220: MatDestroy(&A);
1222: PetscFree(direction);
1223: PetscFree(terminate);
1224: SNESDestroy(&snes_alg);
1225: VecDestroy(&F_alg);
1226: MatDestroy(&J);
1227: MatDestroy(&user.Ybus);
1228: MatDestroy(&user.Sol);
1229: VecDestroy(&X);
1230: VecDestroy(&user.V0);
1231: DMDestroy(&user.dmgen);
1232: DMDestroy(&user.dmnet);
1233: DMDestroy(&user.dmpgrid);
1234: ISDestroy(&user.is_diff);
1235: ISDestroy(&user.is_alg);
1236: TSDestroy(&ts);
1237: if (user.setisdiff) {
1238: VecDestroy(&vatol);
1239: }
1240: PetscFinalize();
1241: return ierr;
1242: }
1244: /*TEST
1246: build:
1247: requires: double !complex !define(PETSC_USE_64BIT_INDICES)
1249: test:
1250: suffix: implicit
1251: args: -ts_monitor -snes_monitor_short
1252: localrunfiles: petscoptions X.bin Ybus.bin
1254: test:
1255: suffix: semiexplicit
1256: args: -ts_monitor -snes_monitor_short -dae_semiexplicit -ts_rk_type 2a
1257: localrunfiles: petscoptions X.bin Ybus.bin
1259: test:
1260: suffix: steprestart
1261: # needs ARKIMEX methods with all implicit stages since the mass matrix is not the identity
1262: args: -ts_monitor -snes_monitor_short -ts_type arkimex -ts_arkimex_type prssp2
1263: localrunfiles: petscoptions X.bin Ybus.bin
1265: TEST*/