Actual source code: ex9busopt_fd.c
petsc-3.6.4 2016-04-12
1: static char help[] = "Power grid stability analysis of WECC 9 bus system.\n\
2: This example is based on the 9-bus (node) example given in the book Power\n\
3: Systems Dynamics and Stability (Chapter 7) by P. Sauer and M. A. Pai.\n\
4: The power grid in this example consists of 9 buses (nodes), 3 generators,\n\
5: 3 loads, and 9 transmission lines. The network equations are written\n\
6: in current balance form using rectangular coordiantes.\n\n";
8: /*
9: The equations for the stability analysis are described by the DAE
11: \dot{x} = f(x,y,t)
12: 0 = g(x,y,t)
14: where the generators are described by differential equations, while the algebraic
15: constraints define the network equations.
17: The generators are modeled with a 4th order differential equation describing the electrical
18: and mechanical dynamics. Each generator also has an exciter system modeled by 3rd order
19: diff. eqns. describing the exciter, voltage regulator, and the feedback stabilizer
20: mechanism.
22: The network equations are described by nodal current balance equations.
23: I(x,y) - Y*V = 0
25: where:
26: I(x,y) is the current injected from generators and loads.
27: Y is the admittance matrix, and
28: V is the voltage vector
29: */
31: /*
32: Include "petscts.h" so that we can use TS solvers. Note that this
33: file automatically includes:
34: petscsys.h - base PETSc routines petscvec.h - vectors
35: petscmat.h - matrices
36: petscis.h - index sets petscksp.h - Krylov subspace methods
37: petscviewer.h - viewers petscpc.h - preconditioners
38: petscksp.h - linear solvers
39: */
40: #include <petsctao.h>
41: #include <petscts.h>
42: #include <petscdm.h>
43: #include <petscdmda.h>
44: #include <petscdmcomposite.h>
46: PetscErrorCode FormFunction(Tao,Vec,PetscReal*,void*);
48: #define freq 60
49: #define w_s (2*PETSC_PI*freq)
51: /* Sizes and indices */
52: const PetscInt nbus = 9; /* Number of network buses */
53: const PetscInt ngen = 3; /* Number of generators */
54: const PetscInt nload = 3; /* Number of loads */
55: const PetscInt gbus[3] = {0,1,2}; /* Buses at which generators are incident */
56: const PetscInt lbus[3] = {4,5,7}; /* Buses at which loads are incident */
58: /* Generator real and reactive powers (found via loadflow) */
59: PetscScalar PG[3] = { 0.69,1.59,0.69};
60: /* PetscScalar PG[3] = {0.716786142395021,1.630000000000000,0.850000000000000};*/
61: const PetscScalar QG[3] = {0.270702180178785,0.066120127797275,-0.108402221791588};
62: /* Generator constants */
63: const PetscScalar H[3] = {23.64,6.4,3.01}; /* Inertia constant */
64: const PetscScalar Rs[3] = {0.0,0.0,0.0}; /* Stator Resistance */
65: const PetscScalar Xd[3] = {0.146,0.8958,1.3125}; /* d-axis reactance */
66: const PetscScalar Xdp[3] = {0.0608,0.1198,0.1813}; /* d-axis transient reactance */
67: 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 */
68: const PetscScalar Xqp[3] = {0.0969,0.1969,0.25}; /* q-axis transient reactance */
69: const PetscScalar Td0p[3] = {8.96,6.0,5.89}; /* d-axis open circuit time constant */
70: const PetscScalar Tq0p[3] = {0.31,0.535,0.6}; /* q-axis open circuit time constant */
71: PetscScalar M[3]; /* M = 2*H/w_s */
72: PetscScalar D[3]; /* D = 0.1*M */
74: PetscScalar TM[3]; /* Mechanical Torque */
75: /* Exciter system constants */
76: const PetscScalar KA[3] = {20.0,20.0,20.0}; /* Voltage regulartor gain constant */
77: const PetscScalar TA[3] = {0.2,0.2,0.2}; /* Voltage regulator time constant */
78: const PetscScalar KE[3] = {1.0,1.0,1.0}; /* Exciter gain constant */
79: const PetscScalar TE[3] = {0.314,0.314,0.314}; /* Exciter time constant */
80: const PetscScalar KF[3] = {0.063,0.063,0.063}; /* Feedback stabilizer gain constant */
81: const PetscScalar TF[3] = {0.35,0.35,0.35}; /* Feedback stabilizer time constant */
82: const PetscScalar k1[3] = {0.0039,0.0039,0.0039};
83: const PetscScalar k2[3] = {1.555,1.555,1.555}; /* k1 and k2 for calculating the saturation function SE = k1*exp(k2*Efd) */
85: PetscScalar Vref[3];
86: /* Load constants
87: We use a composite load model that describes the load and reactive powers at each time instant as follows
88: P(t) = \sum\limits_{i=0}^ld_nsegsp \ld_alphap_i*P_D0(\frac{V_m(t)}{V_m0})^\ld_betap_i
89: Q(t) = \sum\limits_{i=0}^ld_nsegsq \ld_alphaq_i*Q_D0(\frac{V_m(t)}{V_m0})^\ld_betaq_i
90: where
91: ld_nsegsp,ld_nsegsq - Number of individual load models for real and reactive power loads
92: ld_alphap,ld_alphap - Percentage contribution (weights) or loads
93: P_D0 - Real power load
94: Q_D0 - Reactive power load
95: V_m(t) - Voltage magnitude at time t
96: V_m0 - Voltage magnitude at t = 0
97: ld_betap, ld_betaq - exponents describing the load model for real and reactive part
99: Note: All loads have the same characteristic currently.
100: */
101: const PetscScalar PD0[3] = {1.25,0.9,1.0};
102: const PetscScalar QD0[3] = {0.5,0.3,0.35};
103: const PetscInt ld_nsegsp[3] = {3,3,3};
104: const PetscScalar ld_alphap[3] = {1.0,0.0,0.0};
105: const PetscScalar ld_betap[3] = {2.0,1.0,0.0};
106: const PetscInt ld_nsegsq[3] = {3,3,3};
107: const PetscScalar ld_alphaq[3] = {1.0,0.0,0.0};
108: const PetscScalar ld_betaq[3] = {2.0,1.0,0.0};
110: typedef struct {
111: DM dmgen, dmnet; /* DMs to manage generator and network subsystem */
112: DM dmpgrid; /* Composite DM to manage the entire power grid */
113: Mat Ybus; /* Network admittance matrix */
114: Vec V0; /* Initial voltage vector (Power flow solution) */
115: PetscReal tfaulton,tfaultoff; /* Fault on and off times */
116: PetscInt faultbus; /* Fault bus */
117: PetscScalar Rfault;
118: PetscReal t0,tmax;
119: PetscInt neqs_gen,neqs_net,neqs_pgrid;
120: Mat Sol; /* Matrix to save solution at each time step */
121: PetscInt stepnum;
122: PetscBool alg_flg;
123: PetscReal t;
124: IS is_diff; /* indices for differential equations */
125: IS is_alg; /* indices for algebraic equations */
126: PetscReal freq_u,freq_l; /* upper and lower frequency limit */
127: PetscInt pow; /* power coefficient used in the cost function */
128: Vec vec_q;
129: } Userctx;
132: /* Converts from machine frame (dq) to network (phase a real,imag) reference frame */
135: PetscErrorCode dq2ri(PetscScalar Fd,PetscScalar Fq,PetscScalar delta,PetscScalar *Fr, PetscScalar *Fi)
136: {
138: *Fr = Fd*PetscSinScalar(delta) + Fq*PetscCosScalar(delta);
139: *Fi = -Fd*PetscCosScalar(delta) + Fq*PetscSinScalar(delta);
140: return(0);
141: }
143: /* Converts from network frame ([phase a real,imag) to machine (dq) reference frame */
146: PetscErrorCode ri2dq(PetscScalar Fr,PetscScalar Fi,PetscScalar delta,PetscScalar *Fd, PetscScalar *Fq)
147: {
149: *Fd = Fr*PetscSinScalar(delta) - Fi*PetscCosScalar(delta);
150: *Fq = Fr*PetscCosScalar(delta) + Fi*PetscSinScalar(delta);
151: return(0);
152: }
156: PetscErrorCode SetInitialGuess(Vec X,Userctx *user)
157: {
159: Vec Xgen,Xnet;
160: PetscScalar *xgen,*xnet;
161: PetscInt i,idx=0;
162: PetscScalar Vr,Vi,IGr,IGi,Vm,Vm2;
163: PetscScalar Eqp,Edp,delta;
164: PetscScalar Efd,RF,VR; /* Exciter variables */
165: PetscScalar Id,Iq; /* Generator dq axis currents */
166: PetscScalar theta,Vd,Vq,SE;
169: M[0] = 2*H[0]/w_s; M[1] = 2*H[1]/w_s; M[2] = 2*H[2]/w_s;
170: D[0] = 0.1*M[0]; D[1] = 0.1*M[1]; D[2] = 0.1*M[2];
172: DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
174: /* Network subsystem initialization */
175: VecCopy(user->V0,Xnet);
177: /* Generator subsystem initialization */
178: VecGetArray(Xgen,&xgen);
179: VecGetArray(Xnet,&xnet);
181: for (i=0; i < ngen; i++) {
182: Vr = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */
183: Vi = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */
184: Vm = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm;
185: IGr = (Vr*PG[i] + Vi*QG[i])/Vm2;
186: IGi = (Vi*PG[i] - Vr*QG[i])/Vm2;
188: delta = atan2(Vi+Xq[i]*IGr,Vr-Xq[i]*IGi); /* Machine angle */
190: theta = PETSC_PI/2.0 - delta;
192: Id = IGr*PetscCosScalar(theta) - IGi*PetscSinScalar(theta); /* d-axis stator current */
193: Iq = IGr*PetscSinScalar(theta) + IGi*PetscCosScalar(theta); /* q-axis stator current */
195: Vd = Vr*PetscCosScalar(theta) - Vi*PetscSinScalar(theta);
196: Vq = Vr*PetscSinScalar(theta) + Vi*PetscCosScalar(theta);
198: Edp = Vd + Rs[i]*Id - Xqp[i]*Iq; /* d-axis transient EMF */
199: Eqp = Vq + Rs[i]*Iq + Xdp[i]*Id; /* q-axis transient EMF */
201: TM[i] = PG[i];
203: /* The generator variables are ordered as [Eqp,Edp,delta,w,Id,Iq] */
204: xgen[idx] = Eqp;
205: xgen[idx+1] = Edp;
206: xgen[idx+2] = delta;
207: xgen[idx+3] = w_s;
209: idx = idx + 4;
211: xgen[idx] = Id;
212: xgen[idx+1] = Iq;
214: idx = idx + 2;
216: /* Exciter */
217: Efd = Eqp + (Xd[i] - Xdp[i])*Id;
218: SE = k1[i]*PetscExpScalar(k2[i]*Efd);
219: VR = KE[i]*Efd + SE;
220: RF = KF[i]*Efd/TF[i];
222: xgen[idx] = Efd;
223: xgen[idx+1] = RF;
224: xgen[idx+2] = VR;
226: Vref[i] = Vm + (VR/KA[i]);
228: idx = idx + 3;
229: }
231: VecRestoreArray(Xgen,&xgen);
232: VecRestoreArray(Xnet,&xnet);
234: /* VecView(Xgen,0); */
235: DMCompositeGather(user->dmpgrid,X,INSERT_VALUES,Xgen,Xnet);
236: DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);
237: return(0);
238: }
240: /* Computes F = [-f(x,y);g(x,y)] */
243: PetscErrorCode ResidualFunction(SNES snes,Vec X, Vec F, Userctx *user)
244: {
246: Vec Xgen,Xnet,Fgen,Fnet;
247: PetscScalar *xgen,*xnet,*fgen,*fnet;
248: PetscInt i,idx=0;
249: PetscScalar Vr,Vi,Vm,Vm2;
250: PetscScalar Eqp,Edp,delta,w; /* Generator variables */
251: PetscScalar Efd,RF,VR; /* Exciter variables */
252: PetscScalar Id,Iq; /* Generator dq axis currents */
253: PetscScalar Vd,Vq,SE;
254: PetscScalar IGr,IGi,IDr,IDi;
255: PetscScalar Zdq_inv[4],det;
256: PetscScalar PD,QD,Vm0,*v0;
257: PetscInt k;
260: VecZeroEntries(F);
261: DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
262: DMCompositeGetLocalVectors(user->dmpgrid,&Fgen,&Fnet);
263: DMCompositeScatter(user->dmpgrid,X,Xgen,Xnet);
264: DMCompositeScatter(user->dmpgrid,F,Fgen,Fnet);
266: /* Network current balance residual IG + Y*V + IL = 0. Only YV is added here.
267: The generator current injection, IG, and load current injection, ID are added later
268: */
269: /* Note that the values in Ybus are stored assuming the imaginary current balance
270: equation is ordered first followed by real current balance equation for each bus.
271: Thus imaginary current contribution goes in location 2*i, and
272: real current contribution in 2*i+1
273: */
274: MatMult(user->Ybus,Xnet,Fnet);
276: VecGetArray(Xgen,&xgen);
277: VecGetArray(Xnet,&xnet);
278: VecGetArray(Fgen,&fgen);
279: VecGetArray(Fnet,&fnet);
281: /* Generator subsystem */
282: for (i=0; i < ngen; i++) {
283: Eqp = xgen[idx];
284: Edp = xgen[idx+1];
285: delta = xgen[idx+2];
286: w = xgen[idx+3];
287: Id = xgen[idx+4];
288: Iq = xgen[idx+5];
289: Efd = xgen[idx+6];
290: RF = xgen[idx+7];
291: VR = xgen[idx+8];
293: /* Generator differential equations */
294: fgen[idx] = (Eqp + (Xd[i] - Xdp[i])*Id - Efd)/Td0p[i];
295: fgen[idx+1] = (Edp - (Xq[i] - Xqp[i])*Iq)/Tq0p[i];
296: fgen[idx+2] = -w + w_s;
297: fgen[idx+3] = (-TM[i] + Edp*Id + Eqp*Iq + (Xqp[i] - Xdp[i])*Id*Iq + D[i]*(w - w_s))/M[i];
299: Vr = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */
300: Vi = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */
302: ri2dq(Vr,Vi,delta,&Vd,&Vq);
303: /* Algebraic equations for stator currents */
304: det = Rs[i]*Rs[i] + Xdp[i]*Xqp[i];
306: Zdq_inv[0] = Rs[i]/det;
307: Zdq_inv[1] = Xqp[i]/det;
308: Zdq_inv[2] = -Xdp[i]/det;
309: Zdq_inv[3] = Rs[i]/det;
311: fgen[idx+4] = Zdq_inv[0]*(-Edp + Vd) + Zdq_inv[1]*(-Eqp + Vq) + Id;
312: fgen[idx+5] = Zdq_inv[2]*(-Edp + Vd) + Zdq_inv[3]*(-Eqp + Vq) + Iq;
314: /* Add generator current injection to network */
315: dq2ri(Id,Iq,delta,&IGr,&IGi);
317: fnet[2*gbus[i]] -= IGi;
318: fnet[2*gbus[i]+1] -= IGr;
320: Vm = PetscSqrtScalar(Vd*Vd + Vq*Vq); Vm2 = Vm*Vm;
322: SE = k1[i]*PetscExpScalar(k2[i]*Efd);
324: /* Exciter differential equations */
325: fgen[idx+6] = (KE[i]*Efd + SE - VR)/TE[i];
326: fgen[idx+7] = (RF - KF[i]*Efd/TF[i])/TF[i];
327: fgen[idx+8] = (VR - KA[i]*RF + KA[i]*KF[i]*Efd/TF[i] - KA[i]*(Vref[i] - Vm))/TA[i];
329: idx = idx + 9;
330: }
332: VecGetArray(user->V0,&v0);
333: for (i=0; i < nload; i++) {
334: Vr = xnet[2*lbus[i]]; /* Real part of load bus voltage */
335: Vi = xnet[2*lbus[i]+1]; /* Imaginary part of the load bus voltage */
336: Vm = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm;
337: Vm0 = PetscSqrtScalar(v0[2*lbus[i]]*v0[2*lbus[i]] + v0[2*lbus[i]+1]*v0[2*lbus[i]+1]);
338: PD = QD = 0.0;
339: for (k=0; k < ld_nsegsp[i]; k++) PD += ld_alphap[k]*PD0[i]*PetscPowScalar((Vm/Vm0),ld_betap[k]);
340: for (k=0; k < ld_nsegsq[i]; k++) QD += ld_alphaq[k]*QD0[i]*PetscPowScalar((Vm/Vm0),ld_betaq[k]);
342: /* Load currents */
343: IDr = (PD*Vr + QD*Vi)/Vm2;
344: IDi = (-QD*Vr + PD*Vi)/Vm2;
346: fnet[2*lbus[i]] += IDi;
347: fnet[2*lbus[i]+1] += IDr;
348: }
349: VecRestoreArray(user->V0,&v0);
351: VecRestoreArray(Xgen,&xgen);
352: VecRestoreArray(Xnet,&xnet);
353: VecRestoreArray(Fgen,&fgen);
354: VecRestoreArray(Fnet,&fnet);
356: DMCompositeGather(user->dmpgrid,F,INSERT_VALUES,Fgen,Fnet);
357: DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);
358: DMCompositeRestoreLocalVectors(user->dmpgrid,&Fgen,&Fnet);
359: return(0);
360: }
362: /* \dot{x} - f(x,y)
363: g(x,y) = 0
364: */
367: PetscErrorCode IFunction(TS ts,PetscReal t, Vec X, Vec Xdot, Vec F, Userctx *user)
368: {
369: PetscErrorCode ierr;
370: SNES snes;
371: PetscScalar *f;
372: const PetscScalar *xdot;
373: PetscInt i;
376: user->t = t;
378: TSGetSNES(ts,&snes);
379: ResidualFunction(snes,X,F,user);
380: VecGetArray(F,&f);
381: VecGetArrayRead(Xdot,&xdot);
382: for (i=0;i < ngen;i++) {
383: f[9*i] += xdot[9*i];
384: f[9*i+1] += xdot[9*i+1];
385: f[9*i+2] += xdot[9*i+2];
386: f[9*i+3] += xdot[9*i+3];
387: f[9*i+6] += xdot[9*i+6];
388: f[9*i+7] += xdot[9*i+7];
389: f[9*i+8] += xdot[9*i+8];
390: }
391: VecRestoreArray(F,&f);
392: VecRestoreArrayRead(Xdot,&xdot);
393: return(0);
394: }
396: /* This function is used for solving the algebraic system only during fault on and
397: off times. It computes the entire F and then zeros out the part corresponding to
398: differential equations
399: F = [0;g(y)];
400: */
403: PetscErrorCode AlgFunction(SNES snes, Vec X, Vec F, void *ctx)
404: {
406: Userctx *user=(Userctx*)ctx;
407: PetscScalar *f;
408: PetscInt i;
411: ResidualFunction(snes,X,F,user);
412: VecGetArray(F,&f);
413: for (i=0; i < ngen; i++) {
414: f[9*i] = 0;
415: f[9*i+1] = 0;
416: f[9*i+2] = 0;
417: f[9*i+3] = 0;
418: f[9*i+6] = 0;
419: f[9*i+7] = 0;
420: f[9*i+8] = 0;
421: }
422: VecRestoreArray(F,&f);
423: return(0);
424: }
428: PetscErrorCode PreallocateJacobian(Mat J, Userctx *user)
429: {
431: PetscInt *d_nnz;
432: PetscInt i,idx=0,start=0;
433: PetscInt ncols;
436: PetscMalloc1(user->neqs_pgrid,&d_nnz);
437: for (i=0; i<user->neqs_pgrid; i++) d_nnz[i] = 0;
438: /* Generator subsystem */
439: for (i=0; i < ngen; i++) {
441: d_nnz[idx] += 3;
442: d_nnz[idx+1] += 2;
443: d_nnz[idx+2] += 2;
444: d_nnz[idx+3] += 5;
445: d_nnz[idx+4] += 6;
446: d_nnz[idx+5] += 6;
448: d_nnz[user->neqs_gen+2*gbus[i]] += 3;
449: d_nnz[user->neqs_gen+2*gbus[i]+1] += 3;
451: d_nnz[idx+6] += 2;
452: d_nnz[idx+7] += 2;
453: d_nnz[idx+8] += 5;
455: idx = idx + 9;
456: }
458: start = user->neqs_gen;
460: for (i=0; i < nbus; i++) {
461: MatGetRow(user->Ybus,2*i,&ncols,NULL,NULL);
462: d_nnz[start+2*i] += ncols;
463: d_nnz[start+2*i+1] += ncols;
464: MatRestoreRow(user->Ybus,2*i,&ncols,NULL,NULL);
465: }
467: MatSeqAIJSetPreallocation(J,0,d_nnz);
469: PetscFree(d_nnz);
470: return(0);
471: }
473: /*
474: J = [-df_dx, -df_dy
475: dg_dx, dg_dy]
476: */
479: PetscErrorCode ResidualJacobian(SNES snes,Vec X,Mat J,Mat B,void *ctx)
480: {
481: PetscErrorCode ierr;
482: Userctx *user=(Userctx*)ctx;
483: Vec Xgen,Xnet;
484: PetscScalar *xgen,*xnet;
485: PetscInt i,idx=0;
486: PetscScalar Vr,Vi,Vm,Vm2;
487: PetscScalar Eqp,Edp,delta; /* Generator variables */
488: PetscScalar Efd; /* Exciter variables */
489: PetscScalar Id,Iq; /* Generator dq axis currents */
490: PetscScalar Vd,Vq;
491: PetscScalar val[10];
492: PetscInt row[2],col[10];
493: PetscInt net_start=user->neqs_gen;
494: PetscScalar Zdq_inv[4],det;
495: PetscScalar dVd_dVr,dVd_dVi,dVq_dVr,dVq_dVi,dVd_ddelta,dVq_ddelta;
496: PetscScalar dIGr_ddelta,dIGi_ddelta,dIGr_dId,dIGr_dIq,dIGi_dId,dIGi_dIq;
497: PetscScalar dSE_dEfd;
498: PetscScalar dVm_dVd,dVm_dVq,dVm_dVr,dVm_dVi;
499: PetscInt ncols;
500: const PetscInt *cols;
501: const PetscScalar *yvals;
502: PetscInt k;
503: PetscScalar PD,QD,Vm0,*v0,Vm4;
504: PetscScalar dPD_dVr,dPD_dVi,dQD_dVr,dQD_dVi;
505: PetscScalar dIDr_dVr,dIDr_dVi,dIDi_dVr,dIDi_dVi;
508: MatZeroEntries(B);
509: DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
510: DMCompositeScatter(user->dmpgrid,X,Xgen,Xnet);
512: VecGetArray(Xgen,&xgen);
513: VecGetArray(Xnet,&xnet);
515: /* Generator subsystem */
516: for (i=0; i < ngen; i++) {
517: Eqp = xgen[idx];
518: Edp = xgen[idx+1];
519: delta = xgen[idx+2];
520: Id = xgen[idx+4];
521: Iq = xgen[idx+5];
522: Efd = xgen[idx+6];
524: /* fgen[idx] = (Eqp + (Xd[i] - Xdp[i])*Id - Efd)/Td0p[i]; */
525: row[0] = idx;
526: col[0] = idx; col[1] = idx+4; col[2] = idx+6;
527: val[0] = 1/ Td0p[i]; val[1] = (Xd[i] - Xdp[i])/ Td0p[i]; val[2] = -1/Td0p[i];
529: MatSetValues(J,1,row,3,col,val,INSERT_VALUES);
531: /* fgen[idx+1] = (Edp - (Xq[i] - Xqp[i])*Iq)/Tq0p[i]; */
532: row[0] = idx + 1;
533: col[0] = idx + 1; col[1] = idx+5;
534: val[0] = 1/Tq0p[i]; val[1] = -(Xq[i] - Xqp[i])/Tq0p[i];
535: MatSetValues(J,1,row,2,col,val,INSERT_VALUES);
537: /* fgen[idx+2] = - w + w_s; */
538: row[0] = idx + 2;
539: col[0] = idx + 2; col[1] = idx + 3;
540: val[0] = 0; val[1] = -1;
541: MatSetValues(J,1,row,2,col,val,INSERT_VALUES);
543: /* fgen[idx+3] = (-TM[i] + Edp*Id + Eqp*Iq + (Xqp[i] - Xdp[i])*Id*Iq + D[i]*(w - w_s))/M[i]; */
544: row[0] = idx + 3;
545: col[0] = idx; col[1] = idx + 1; col[2] = idx + 3; col[3] = idx + 4; col[4] = idx + 5;
546: 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];
547: MatSetValues(J,1,row,5,col,val,INSERT_VALUES);
549: Vr = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */
550: Vi = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */
551: ri2dq(Vr,Vi,delta,&Vd,&Vq);
553: det = Rs[i]*Rs[i] + Xdp[i]*Xqp[i];
555: Zdq_inv[0] = Rs[i]/det;
556: Zdq_inv[1] = Xqp[i]/det;
557: Zdq_inv[2] = -Xdp[i]/det;
558: Zdq_inv[3] = Rs[i]/det;
560: dVd_dVr = PetscSinScalar(delta); dVd_dVi = -PetscCosScalar(delta);
561: dVq_dVr = PetscCosScalar(delta); dVq_dVi = PetscSinScalar(delta);
562: dVd_ddelta = Vr*PetscCosScalar(delta) + Vi*PetscSinScalar(delta);
563: dVq_ddelta = -Vr*PetscSinScalar(delta) + Vi*PetscCosScalar(delta);
565: /* fgen[idx+4] = Zdq_inv[0]*(-Edp + Vd) + Zdq_inv[1]*(-Eqp + Vq) + Id; */
566: row[0] = idx+4;
567: col[0] = idx; col[1] = idx+1; col[2] = idx + 2;
568: val[0] = -Zdq_inv[1]; val[1] = -Zdq_inv[0]; val[2] = Zdq_inv[0]*dVd_ddelta + Zdq_inv[1]*dVq_ddelta;
569: col[3] = idx + 4; col[4] = net_start+2*gbus[i]; col[5] = net_start + 2*gbus[i]+1;
570: 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;
571: MatSetValues(J,1,row,6,col,val,INSERT_VALUES);
573: /* fgen[idx+5] = Zdq_inv[2]*(-Edp + Vd) + Zdq_inv[3]*(-Eqp + Vq) + Iq; */
574: row[0] = idx+5;
575: col[0] = idx; col[1] = idx+1; col[2] = idx + 2;
576: val[0] = -Zdq_inv[3]; val[1] = -Zdq_inv[2]; val[2] = Zdq_inv[2]*dVd_ddelta + Zdq_inv[3]*dVq_ddelta;
577: col[3] = idx + 5; col[4] = net_start+2*gbus[i]; col[5] = net_start + 2*gbus[i]+1;
578: 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;
579: MatSetValues(J,1,row,6,col,val,INSERT_VALUES);
581: dIGr_ddelta = Id*PetscCosScalar(delta) - Iq*PetscSinScalar(delta);
582: dIGi_ddelta = Id*PetscSinScalar(delta) + Iq*PetscCosScalar(delta);
583: dIGr_dId = PetscSinScalar(delta); dIGr_dIq = PetscCosScalar(delta);
584: dIGi_dId = -PetscCosScalar(delta); dIGi_dIq = PetscSinScalar(delta);
586: /* fnet[2*gbus[i]] -= IGi; */
587: row[0] = net_start + 2*gbus[i];
588: col[0] = idx+2; col[1] = idx + 4; col[2] = idx + 5;
589: val[0] = -dIGi_ddelta; val[1] = -dIGi_dId; val[2] = -dIGi_dIq;
590: MatSetValues(J,1,row,3,col,val,INSERT_VALUES);
592: /* fnet[2*gbus[i]+1] -= IGr; */
593: row[0] = net_start + 2*gbus[i]+1;
594: col[0] = idx+2; col[1] = idx + 4; col[2] = idx + 5;
595: val[0] = -dIGr_ddelta; val[1] = -dIGr_dId; val[2] = -dIGr_dIq;
596: MatSetValues(J,1,row,3,col,val,INSERT_VALUES);
598: Vm = PetscSqrtScalar(Vd*Vd + Vq*Vq); Vm2 = Vm*Vm;
600: /* fgen[idx+6] = (KE[i]*Efd + SE - VR)/TE[i]; */
601: /* SE = k1[i]*PetscExpScalar(k2[i]*Efd); */
603: dSE_dEfd = k1[i]*k2[i]*PetscExpScalar(k2[i]*Efd);
605: row[0] = idx + 6;
606: col[0] = idx + 6; col[1] = idx + 8;
607: val[0] = (KE[i] + dSE_dEfd)/TE[i]; val[1] = -1/TE[i];
608: MatSetValues(J,1,row,2,col,val,INSERT_VALUES);
610: /* Exciter differential equations */
612: /* fgen[idx+7] = (RF - KF[i]*Efd/TF[i])/TF[i]; */
613: row[0] = idx + 7;
614: col[0] = idx + 6; col[1] = idx + 7;
615: val[0] = (-KF[i]/TF[i])/TF[i]; val[1] = 1/TF[i];
616: MatSetValues(J,1,row,2,col,val,INSERT_VALUES);
618: /* fgen[idx+8] = (VR - KA[i]*RF + KA[i]*KF[i]*Efd/TF[i] - KA[i]*(Vref[i] - Vm))/TA[i]; */
619: /* Vm = (Vd^2 + Vq^2)^0.5; */
620: dVm_dVd = Vd/Vm; dVm_dVq = Vq/Vm;
621: dVm_dVr = dVm_dVd*dVd_dVr + dVm_dVq*dVq_dVr;
622: dVm_dVi = dVm_dVd*dVd_dVi + dVm_dVq*dVq_dVi;
623: row[0] = idx + 8;
624: col[0] = idx + 6; col[1] = idx + 7; col[2] = idx + 8;
625: val[0] = (KA[i]*KF[i]/TF[i])/TA[i]; val[1] = -KA[i]/TA[i]; val[2] = 1/TA[i];
626: col[3] = net_start + 2*gbus[i]; col[4] = net_start + 2*gbus[i]+1;
627: val[3] = KA[i]*dVm_dVr/TA[i]; val[4] = KA[i]*dVm_dVi/TA[i];
628: MatSetValues(J,1,row,5,col,val,INSERT_VALUES);
629: idx = idx + 9;
630: }
633: for (i=0; i<nbus; i++) {
634: MatGetRow(user->Ybus,2*i,&ncols,&cols,&yvals);
635: row[0] = net_start + 2*i;
636: for (k=0; k<ncols; k++) {
637: col[k] = net_start + cols[k];
638: val[k] = yvals[k];
639: }
640: MatSetValues(J,1,row,ncols,col,val,INSERT_VALUES);
641: MatRestoreRow(user->Ybus,2*i,&ncols,&cols,&yvals);
643: MatGetRow(user->Ybus,2*i+1,&ncols,&cols,&yvals);
644: row[0] = net_start + 2*i+1;
645: for (k=0; k<ncols; k++) {
646: col[k] = net_start + cols[k];
647: val[k] = yvals[k];
648: }
649: MatSetValues(J,1,row,ncols,col,val,INSERT_VALUES);
650: MatRestoreRow(user->Ybus,2*i+1,&ncols,&cols,&yvals);
651: }
653: MatAssemblyBegin(J,MAT_FLUSH_ASSEMBLY);
654: MatAssemblyEnd(J,MAT_FLUSH_ASSEMBLY);
656: VecGetArray(user->V0,&v0);
657: for (i=0; i < nload; i++) {
658: Vr = xnet[2*lbus[i]]; /* Real part of load bus voltage */
659: Vi = xnet[2*lbus[i]+1]; /* Imaginary part of the load bus voltage */
660: Vm = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm; Vm4 = Vm2*Vm2;
661: Vm0 = PetscSqrtScalar(v0[2*lbus[i]]*v0[2*lbus[i]] + v0[2*lbus[i]+1]*v0[2*lbus[i]+1]);
662: PD = QD = 0.0;
663: dPD_dVr = dPD_dVi = dQD_dVr = dQD_dVi = 0.0;
664: for (k=0; k < ld_nsegsp[i]; k++) {
665: PD += ld_alphap[k]*PD0[i]*PetscPowScalar((Vm/Vm0),ld_betap[k]);
666: dPD_dVr += ld_alphap[k]*ld_betap[k]*PD0[i]*PetscPowScalar((1/Vm0),ld_betap[k])*Vr*PetscPowScalar(Vm,(ld_betap[k]-2));
667: dPD_dVi += ld_alphap[k]*ld_betap[k]*PD0[i]*PetscPowScalar((1/Vm0),ld_betap[k])*Vi*PetscPowScalar(Vm,(ld_betap[k]-2));
668: }
669: for (k=0; k < ld_nsegsq[i]; k++) {
670: QD += ld_alphaq[k]*QD0[i]*PetscPowScalar((Vm/Vm0),ld_betaq[k]);
671: dQD_dVr += ld_alphaq[k]*ld_betaq[k]*QD0[i]*PetscPowScalar((1/Vm0),ld_betaq[k])*Vr*PetscPowScalar(Vm,(ld_betaq[k]-2));
672: dQD_dVi += ld_alphaq[k]*ld_betaq[k]*QD0[i]*PetscPowScalar((1/Vm0),ld_betaq[k])*Vi*PetscPowScalar(Vm,(ld_betaq[k]-2));
673: }
675: /* IDr = (PD*Vr + QD*Vi)/Vm2; */
676: /* IDi = (-QD*Vr + PD*Vi)/Vm2; */
678: dIDr_dVr = (dPD_dVr*Vr + dQD_dVr*Vi + PD)/Vm2 - ((PD*Vr + QD*Vi)*2*Vr)/Vm4;
679: dIDr_dVi = (dPD_dVi*Vr + dQD_dVi*Vi + QD)/Vm2 - ((PD*Vr + QD*Vi)*2*Vi)/Vm4;
681: dIDi_dVr = (-dQD_dVr*Vr + dPD_dVr*Vi - QD)/Vm2 - ((-QD*Vr + PD*Vi)*2*Vr)/Vm4;
682: dIDi_dVi = (-dQD_dVi*Vr + dPD_dVi*Vi + PD)/Vm2 - ((-QD*Vr + PD*Vi)*2*Vi)/Vm4;
685: /* fnet[2*lbus[i]] += IDi; */
686: row[0] = net_start + 2*lbus[i];
687: col[0] = net_start + 2*lbus[i]; col[1] = net_start + 2*lbus[i]+1;
688: val[0] = dIDi_dVr; val[1] = dIDi_dVi;
689: MatSetValues(J,1,row,2,col,val,ADD_VALUES);
690: /* fnet[2*lbus[i]+1] += IDr; */
691: row[0] = net_start + 2*lbus[i]+1;
692: col[0] = net_start + 2*lbus[i]; col[1] = net_start + 2*lbus[i]+1;
693: val[0] = dIDr_dVr; val[1] = dIDr_dVi;
694: MatSetValues(J,1,row,2,col,val,ADD_VALUES);
695: }
696: VecRestoreArray(user->V0,&v0);
698: VecRestoreArray(Xgen,&xgen);
699: VecRestoreArray(Xnet,&xnet);
701: DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);
703: MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
704: MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
705: return(0);
706: }
708: /*
709: J = [I, 0
710: dg_dx, dg_dy]
711: */
714: PetscErrorCode AlgJacobian(SNES snes,Vec X,Mat A,Mat B,void *ctx)
715: {
717: Userctx *user=(Userctx*)ctx;
720: ResidualJacobian(snes,X,A,B,ctx);
721: MatSetOption(A,MAT_KEEP_NONZERO_PATTERN,PETSC_TRUE);
722: MatZeroRowsIS(A,user->is_diff,1.0,NULL,NULL);
723: return(0);
724: }
726: /*
727: J = [a*I-df_dx, -df_dy
728: dg_dx, dg_dy]
729: */
733: PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,Userctx *user)
734: {
736: SNES snes;
737: PetscScalar atmp = (PetscScalar) a;
738: PetscInt i,row;
741: user->t = t;
743: TSGetSNES(ts,&snes);
744: ResidualJacobian(snes,X,A,B,user);
745: for (i=0;i < ngen;i++) {
746: row = 9*i;
747: MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
748: row = 9*i+1;
749: MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
750: row = 9*i+2;
751: MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
752: row = 9*i+3;
753: MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
754: row = 9*i+6;
755: MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
756: row = 9*i+7;
757: MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
758: row = 9*i+8;
759: MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
760: }
761: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
762: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
763: return(0);
764: }
768: static PetscErrorCode CostIntegrand(TS ts,PetscReal t,Vec U,Vec R,Userctx *user)
769: {
770: PetscErrorCode ierr;
771: PetscScalar *r;
772: const PetscScalar *u;
773: PetscInt idx;
774: Vec Xgen,Xnet;
775: PetscScalar *xgen;
776: PetscInt i;
779: DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
780: DMCompositeScatter(user->dmpgrid,U,Xgen,Xnet);
782: VecGetArray(Xgen,&xgen);
784: VecGetArrayRead(U,&u);
785: VecGetArray(R,&r);
786: r[0] = 0.;
788: idx = 0;
789: for (i=0;i<ngen;i++) {
790: r[0] += PetscPowScalarInt(PetscMax(0.,PetscMax(xgen[idx+3]/(2.*PETSC_PI)-user->freq_u,user->freq_l-xgen[idx+3]/(2.*PETSC_PI))),user->pow);
791: idx += 9;
792: }
793: VecRestoreArray(R,&r);
794: VecRestoreArrayRead(U,&u);
795: DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);
796: return(0);
797: }
801: static PetscErrorCode MonitorUpdateQ(TS ts,PetscInt stepnum,PetscReal time,Vec X,void *ctx0)
802: {
804: Vec C,*Y;
805: PetscInt Nr;
806: PetscReal h,theta;
807: Userctx *ctx=(Userctx*)ctx0;
808:
810: theta = 0.5;
811: TSGetStages(ts,&Nr,&Y);
812: TSGetTimeStep(ts,&h);
813: VecDuplicate(ctx->vec_q,&C);
814: /* compute integrals */
815: if (stepnum>0) {
816: CostIntegrand(ts,time,X,C,ctx);
817: VecAXPY(ctx->vec_q,h*theta,C);
818: CostIntegrand(ts,time+h*theta,Y[0],C,ctx);
819: VecAXPY(ctx->vec_q,h*(1-theta),C);
820: }
821: VecDestroy(&C);
822: return(0);
823: }
827: int main(int argc,char **argv)
828: {
829: Userctx user;
830: Vec p;
831: PetscScalar *x_ptr;
832: PetscErrorCode ierr;
833: PetscMPIInt size;
834: PetscInt i;
835: KSP ksp;
836: PC pc;
837: PetscInt *idx2;
838: Tao tao;
839: TaoConvergedReason reason;
840: Vec lowerb,upperb;
842: PetscInitialize(&argc,&argv,"petscoptions",help);
844: MPI_Comm_size(PETSC_COMM_WORLD,&size);
845: if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs");
847: VecCreateSeq(PETSC_COMM_WORLD,1,&user.vec_q);
849: user.neqs_gen = 9*ngen; /* # eqs. for generator subsystem */
850: user.neqs_net = 2*nbus; /* # eqs. for network subsystem */
851: user.neqs_pgrid = user.neqs_gen + user.neqs_net;
853: /* Create indices for differential and algebraic equations */
854: PetscMalloc1(7*ngen,&idx2);
855: for (i=0; i<ngen; i++) {
856: 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;
857: idx2[7*i+4] = 9*i+6; idx2[7*i+5] = 9*i+7; idx2[7*i+6] = 9*i+8;
858: }
859: ISCreateGeneral(PETSC_COMM_WORLD,7*ngen,idx2,PETSC_COPY_VALUES,&user.is_diff);
860: ISComplement(user.is_diff,0,user.neqs_pgrid,&user.is_alg);
861: PetscFree(idx2);
863: /* Set run time options */
864: PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Transient stability fault options","");
865: {
866: user.tfaulton = 1.0;
867: user.tfaultoff = 1.2;
868: user.Rfault = 0.0001;
869: user.faultbus = 8;
870: PetscOptionsReal("-tfaulton","","",user.tfaulton,&user.tfaulton,NULL);
871: PetscOptionsReal("-tfaultoff","","",user.tfaultoff,&user.tfaultoff,NULL);
872: PetscOptionsInt("-faultbus","","",user.faultbus,&user.faultbus,NULL);
873: user.t0 = 0.0;
874: user.tmax = 5.0;
875: PetscOptionsReal("-t0","","",user.t0,&user.t0,NULL);
876: PetscOptionsReal("-tmax","","",user.tmax,&user.tmax,NULL);
877: user.freq_u = 61.0;
878: user.freq_l = 59.0;
879: user.pow = 2;
880: PetscOptionsReal("-frequ","","",user.freq_u,&user.freq_u,NULL);
881: PetscOptionsReal("-freql","","",user.freq_l,&user.freq_l,NULL);
882: PetscOptionsInt("-pow","","",user.pow,&user.pow,NULL);
884: }
885: PetscOptionsEnd();
887: /* Create DMs for generator and network subsystems */
888: DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,user.neqs_gen,1,1,NULL,&user.dmgen);
889: DMSetOptionsPrefix(user.dmgen,"dmgen_");
890: DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,user.neqs_net,1,1,NULL,&user.dmnet);
891: DMSetOptionsPrefix(user.dmnet,"dmnet_");
892: /* Create a composite DM packer and add the two DMs */
893: DMCompositeCreate(PETSC_COMM_WORLD,&user.dmpgrid);
894: DMSetOptionsPrefix(user.dmpgrid,"pgrid_");
895: DMCompositeAddDM(user.dmpgrid,user.dmgen);
896: DMCompositeAddDM(user.dmpgrid,user.dmnet);
898: /* Create TAO solver and set desired solution method */
899: TaoCreate(PETSC_COMM_WORLD,&tao);
900: TaoSetType(tao,TAOBLMVM);
901: /*
902: Optimization starts
903: */
904: /* Set initial solution guess */
905: VecCreateSeq(PETSC_COMM_WORLD,3,&p);
906: VecGetArray(p,&x_ptr);
907: x_ptr[0] = PG[0]; x_ptr[1] = PG[1]; x_ptr[2] = PG[2];
908: VecRestoreArray(p,&x_ptr);
910: TaoSetInitialVector(tao,p);
911: /* Set routine for function and gradient evaluation */
912: TaoSetObjectiveRoutine(tao,FormFunction,(void *)&user);
913: TaoSetGradientRoutine(tao,TaoDefaultComputeGradient,(void *)&user);
915: /* Set bounds for the optimization */
916: VecDuplicate(p,&lowerb);
917: VecDuplicate(p,&upperb);
918: VecGetArray(lowerb,&x_ptr);
919: x_ptr[0] = 0.5; x_ptr[1] = 0.5; x_ptr[2] = 0.5;
920: VecRestoreArray(lowerb,&x_ptr);
921: VecGetArray(upperb,&x_ptr);
922: x_ptr[0] = 2.0; x_ptr[1] = 2.0; x_ptr[2] = 2.0;
923: VecRestoreArray(upperb,&x_ptr);
924: TaoSetVariableBounds(tao,lowerb,upperb);
926: /* Check for any TAO command line options */
927: TaoSetFromOptions(tao);
928: TaoGetKSP(tao,&ksp);
929: if (ksp) {
930: KSPGetPC(ksp,&pc);
931: PCSetType(pc,PCNONE);
932: }
934: /* TaoSetTolerances(tao,1e-15,1e-15,1e-15,1e-15,1e-15); */
935: /* SOLVE THE APPLICATION */
936: TaoSolve(tao);
937: /* Get information on termination */
938: TaoGetConvergedReason(tao,&reason);
939: if (reason <= 0){
940: ierr=PetscPrintf(MPI_COMM_WORLD, "Try another method! \n");
941: }
943: VecView(p,PETSC_VIEWER_STDOUT_WORLD);
944: /* Free TAO data structures */
945: TaoDestroy(&tao);
947: DMDestroy(&user.dmgen);
948: DMDestroy(&user.dmnet);
949: DMDestroy(&user.dmpgrid);
950: ISDestroy(&user.is_diff);
951: ISDestroy(&user.is_alg);
953: return(0);
954: }
956: /* ------------------------------------------------------------------ */
959: /*
960: FormFunction - Evaluates the function and corresponding gradient.
962: Input Parameters:
963: tao - the Tao context
964: X - the input vector
965: ptr - optional user-defined context, as set by TaoSetObjectiveAndGradientRoutine()
967: Output Parameters:
968: f - the newly evaluated function
969: */
970: PetscErrorCode FormFunction(Tao tao,Vec P,PetscReal *f,void *ctx0)
971: {
972: TS ts;
973: SNES snes_alg;
975: Userctx *ctx = (Userctx*)ctx0;
976: Vec X;
977: Mat J;
978: /* sensitivity context */
979: PetscScalar *x_ptr;
980: PetscViewer Xview,Ybusview;
981: Vec F_alg;
982: Vec Xdot;
983: PetscInt row_loc,col_loc;
984: PetscScalar val;
986: VecGetArray(P,&x_ptr);
987: PG[0] = x_ptr[0];
988: PG[1] = x_ptr[1];
989: PG[2] = x_ptr[2];
990: VecRestoreArray(P,&x_ptr);
992: ctx->stepnum = 0;
994: VecZeroEntries(ctx->vec_q);
996: /* Read initial voltage vector and Ybus */
997: PetscViewerBinaryOpen(PETSC_COMM_WORLD,"X.bin",FILE_MODE_READ,&Xview);
998: PetscViewerBinaryOpen(PETSC_COMM_WORLD,"Ybus.bin",FILE_MODE_READ,&Ybusview);
1000: VecCreate(PETSC_COMM_WORLD,&ctx->V0);
1001: VecSetSizes(ctx->V0,PETSC_DECIDE,ctx->neqs_net);
1002: VecLoad(ctx->V0,Xview);
1004: MatCreate(PETSC_COMM_WORLD,&ctx->Ybus);
1005: MatSetSizes(ctx->Ybus,PETSC_DECIDE,PETSC_DECIDE,ctx->neqs_net,ctx->neqs_net);
1006: MatSetType(ctx->Ybus,MATBAIJ);
1007: /* MatSetBlockSize(ctx->Ybus,2); */
1008: MatLoad(ctx->Ybus,Ybusview);
1010: PetscViewerDestroy(&Xview);
1011: PetscViewerDestroy(&Ybusview);
1013: DMCreateGlobalVector(ctx->dmpgrid,&X);
1015: MatCreate(PETSC_COMM_WORLD,&J);
1016: MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,ctx->neqs_pgrid,ctx->neqs_pgrid);
1017: MatSetFromOptions(J);
1018: PreallocateJacobian(J,ctx);
1020: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1021: Create timestepping solver context
1022: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1023: TSCreate(PETSC_COMM_WORLD,&ts);
1024: TSSetProblemType(ts,TS_NONLINEAR);
1025: TSSetType(ts,TSCN);
1026: TSSetIFunction(ts,NULL,(TSIFunction) IFunction,ctx);
1027: TSSetIJacobian(ts,J,J,(TSIJacobian)IJacobian,ctx);
1028: TSSetApplicationContext(ts,ctx);
1030: TSMonitorSet(ts,MonitorUpdateQ,ctx,NULL);
1031: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1032: Set initial conditions
1033: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1034: SetInitialGuess(X,ctx);
1036: VecDuplicate(X,&F_alg);
1037: SNESCreate(PETSC_COMM_WORLD,&snes_alg);
1038: SNESSetFunction(snes_alg,F_alg,AlgFunction,ctx);
1039: MatZeroEntries(J);
1040: SNESSetJacobian(snes_alg,J,J,AlgJacobian,ctx);
1041: SNESSetOptionsPrefix(snes_alg,"alg_");
1042: SNESSetFromOptions(snes_alg);
1043: ctx->alg_flg = PETSC_TRUE;
1044: /* Solve the algebraic equations */
1045: SNESSolve(snes_alg,NULL,X);
1047: /* Just to set up the Jacobian structure */
1048: VecDuplicate(X,&Xdot);
1049: IJacobian(ts,0.0,X,Xdot,0.0,J,J,ctx);
1050: VecDestroy(&Xdot);
1052: ctx->stepnum++;
1054: TSSetDuration(ts,1000,ctx->tfaulton);
1055: TSSetInitialTimeStep(ts,0.0,0.01);
1056: TSSetFromOptions(ts);
1057: /* TSSetPostStep(ts,SaveSolution); */
1059: ctx->alg_flg = PETSC_FALSE;
1060: /* Prefault period */
1061: TSSolve(ts,X);
1063: /* Create the nonlinear solver for solving the algebraic system */
1064: /* Note that although the algebraic system needs to be solved only for
1065: Idq and V, we reuse the entire system including xgen. The xgen
1066: variables are held constant by setting their residuals to 0 and
1067: putting a 1 on the Jacobian diagonal for xgen rows
1068: */
1069: MatZeroEntries(J);
1071: /* Apply disturbance - resistive fault at ctx->faultbus */
1072: /* This is done by adding shunt conductance to the diagonal location
1073: in the Ybus matrix */
1074: row_loc = 2*ctx->faultbus; col_loc = 2*ctx->faultbus+1; /* Location for G */
1075: val = 1/ctx->Rfault;
1076: MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);
1077: row_loc = 2*ctx->faultbus+1; col_loc = 2*ctx->faultbus; /* Location for G */
1078: val = 1/ctx->Rfault;
1079: MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);
1081: MatAssemblyBegin(ctx->Ybus,MAT_FINAL_ASSEMBLY);
1082: MatAssemblyEnd(ctx->Ybus,MAT_FINAL_ASSEMBLY);
1084: ctx->alg_flg = PETSC_TRUE;
1085: /* Solve the algebraic equations */
1086: SNESSolve(snes_alg,NULL,X);
1088: ctx->stepnum++;
1090: /* Disturbance period */
1091: TSSetDuration(ts,1000,ctx->tfaultoff);
1092: TSSetInitialTimeStep(ts,ctx->tfaulton,.01);
1094: ctx->alg_flg = PETSC_FALSE;
1096: TSSolve(ts,X);
1098: /* Remove the fault */
1099: row_loc = 2*ctx->faultbus; col_loc = 2*ctx->faultbus+1;
1100: val = -1/ctx->Rfault;
1101: MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);
1102: row_loc = 2*ctx->faultbus+1; col_loc = 2*ctx->faultbus;
1103: val = -1/ctx->Rfault;
1104: MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);
1106: MatAssemblyBegin(ctx->Ybus,MAT_FINAL_ASSEMBLY);
1107: MatAssemblyEnd(ctx->Ybus,MAT_FINAL_ASSEMBLY);
1109: MatZeroEntries(J);
1111: ctx->alg_flg = PETSC_TRUE;
1113: /* Solve the algebraic equations */
1114: SNESSolve(snes_alg,NULL,X);
1116: ctx->stepnum++;
1118: /* Post-disturbance period */
1119: TSSetDuration(ts,1000,ctx->tmax);
1120: TSSetInitialTimeStep(ts,ctx->tfaultoff,.01);
1122: ctx->alg_flg = PETSC_TRUE;
1124: TSSolve(ts,X);
1125: VecGetArray(ctx->vec_q,&x_ptr);
1126: *f = x_ptr[0];
1127: VecRestoreArray(ctx->vec_q,&x_ptr);
1129: MatDestroy(&ctx->Ybus);
1130: VecDestroy(&ctx->V0);
1131: SNESDestroy(&snes_alg);
1132: VecDestroy(&F_alg);
1133: MatDestroy(&J);
1134: VecDestroy(&X);
1135: TSDestroy(&ts);
1137: return 0;
1138: }