Actual source code: ex9busopt.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>
45: #include <petsctime.h>
47: PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal*,Vec,void*);
49: #define freq 60
50: #define w_s (2*PETSC_PI*freq)
52: /* Sizes and indices */
53: const PetscInt nbus = 9; /* Number of network buses */
54: const PetscInt ngen = 3; /* Number of generators */
55: const PetscInt nload = 3; /* Number of loads */
56: const PetscInt gbus[3] = {0,1,2}; /* Buses at which generators are incident */
57: const PetscInt lbus[3] = {4,5,7}; /* Buses at which loads are incident */
59: /* Generator real and reactive powers (found via loadflow) */
60: PetscScalar PG[3] = { 0.69,1.59,0.69};
61: /* PetscScalar PG[3] = {0.716786142395021,1.630000000000000,0.850000000000000};*/
63: const PetscScalar QG[3] = {0.270702180178785,0.066120127797275,-0.108402221791588};
64: /* Generator constants */
65: const PetscScalar H[3] = {23.64,6.4,3.01}; /* Inertia constant */
66: const PetscScalar Rs[3] = {0.0,0.0,0.0}; /* Stator Resistance */
67: const PetscScalar Xd[3] = {0.146,0.8958,1.3125}; /* d-axis reactance */
68: const PetscScalar Xdp[3] = {0.0608,0.1198,0.1813}; /* d-axis transient reactance */
69: 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 */
70: const PetscScalar Xqp[3] = {0.0969,0.1969,0.25}; /* q-axis transient reactance */
71: const PetscScalar Td0p[3] = {8.96,6.0,5.89}; /* d-axis open circuit time constant */
72: const PetscScalar Tq0p[3] = {0.31,0.535,0.6}; /* q-axis open circuit time constant */
73: PetscScalar M[3]; /* M = 2*H/w_s */
74: PetscScalar D[3]; /* D = 0.1*M */
76: PetscScalar TM[3]; /* Mechanical Torque */
77: /* Exciter system constants */
78: const PetscScalar KA[3] = {20.0,20.0,20.0}; /* Voltage regulartor gain constant */
79: const PetscScalar TA[3] = {0.2,0.2,0.2}; /* Voltage regulator time constant */
80: const PetscScalar KE[3] = {1.0,1.0,1.0}; /* Exciter gain constant */
81: const PetscScalar TE[3] = {0.314,0.314,0.314}; /* Exciter time constant */
82: const PetscScalar KF[3] = {0.063,0.063,0.063}; /* Feedback stabilizer gain constant */
83: const PetscScalar TF[3] = {0.35,0.35,0.35}; /* Feedback stabilizer time constant */
84: const PetscScalar k1[3] = {0.0039,0.0039,0.0039};
85: const PetscScalar k2[3] = {1.555,1.555,1.555}; /* k1 and k2 for calculating the saturation function SE = k1*exp(k2*Efd) */
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: PetscBool alg_flg;
125: PetscReal t;
126: IS is_diff; /* indices for differential equations */
127: IS is_alg; /* indices for algebraic equations */
128: PetscReal freq_u,freq_l; /* upper and lower frequency limit */
129: PetscInt pow; /* power coefficient used in the cost function */
130: PetscBool jacp_flg;
131: Mat J,Jacp;
132: } Userctx;
135: /* Converts from machine frame (dq) to network (phase a real,imag) reference frame */
138: PetscErrorCode dq2ri(PetscScalar Fd,PetscScalar Fq,PetscScalar delta,PetscScalar *Fr, PetscScalar *Fi)
139: {
141: *Fr = Fd*PetscSinScalar(delta) + Fq*PetscCosScalar(delta);
142: *Fi = -Fd*PetscCosScalar(delta) + Fq*PetscSinScalar(delta);
143: return(0);
144: }
146: /* Converts from network frame ([phase a real,imag) to machine (dq) reference frame */
149: PetscErrorCode ri2dq(PetscScalar Fr,PetscScalar Fi,PetscScalar delta,PetscScalar *Fd, PetscScalar *Fq)
150: {
152: *Fd = Fr*PetscSinScalar(delta) - Fi*PetscCosScalar(delta);
153: *Fq = Fr*PetscCosScalar(delta) + Fi*PetscSinScalar(delta);
154: return(0);
155: }
157: /* Saves the solution at each time to a matrix */
160: PetscErrorCode SaveSolution(TS ts)
161: {
162: PetscErrorCode ierr;
163: Userctx *user;
164: Vec X;
165: PetscScalar *mat;
166: const PetscScalar *x;
167: PetscInt idx;
168: PetscReal t;
171: TSGetApplicationContext(ts,&user);
172: TSGetTime(ts,&t);
173: TSGetSolution(ts,&X);
174: idx = user->stepnum*(user->neqs_pgrid+1);
175: MatDenseGetArray(user->Sol,&mat);
176: VecGetArrayRead(X,&x);
177: mat[idx] = t;
178: PetscMemcpy(mat+idx+1,x,user->neqs_pgrid*sizeof(PetscScalar));
179: MatDenseRestoreArray(user->Sol,&mat);
180: VecRestoreArrayRead(X,&x);
181: user->stepnum++;
182: return(0);
183: }
187: PetscErrorCode SetInitialGuess(Vec X,Userctx *user)
188: {
190: Vec Xgen,Xnet;
191: PetscScalar *xgen,*xnet;
192: PetscInt i,idx=0;
193: PetscScalar Vr,Vi,IGr,IGi,Vm,Vm2;
194: PetscScalar Eqp,Edp,delta;
195: PetscScalar Efd,RF,VR; /* Exciter variables */
196: PetscScalar Id,Iq; /* Generator dq axis currents */
197: PetscScalar theta,Vd,Vq,SE;
200: M[0] = 2*H[0]/w_s; M[1] = 2*H[1]/w_s; M[2] = 2*H[2]/w_s;
201: D[0] = 0.1*M[0]; D[1] = 0.1*M[1]; D[2] = 0.1*M[2];
203: DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
205: /* Network subsystem initialization */
206: VecCopy(user->V0,Xnet);
208: /* Generator subsystem initialization */
209: VecGetArray(Xgen,&xgen);
210: VecGetArray(Xnet,&xnet);
212: for (i=0; i < ngen; i++) {
213: Vr = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */
214: Vi = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */
215: Vm = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm;
216: IGr = (Vr*PG[i] + Vi*QG[i])/Vm2;
217: IGi = (Vi*PG[i] - Vr*QG[i])/Vm2;
219: delta = atan2(Vi+Xq[i]*IGr,Vr-Xq[i]*IGi); /* Machine angle */
221: theta = PETSC_PI/2.0 - delta;
223: Id = IGr*PetscCosScalar(theta) - IGi*PetscSinScalar(theta); /* d-axis stator current */
224: Iq = IGr*PetscSinScalar(theta) + IGi*PetscCosScalar(theta); /* q-axis stator current */
226: Vd = Vr*PetscCosScalar(theta) - Vi*PetscSinScalar(theta);
227: Vq = Vr*PetscSinScalar(theta) + Vi*PetscCosScalar(theta);
229: Edp = Vd + Rs[i]*Id - Xqp[i]*Iq; /* d-axis transient EMF */
230: Eqp = Vq + Rs[i]*Iq + Xdp[i]*Id; /* q-axis transient EMF */
232: TM[i] = PG[i];
234: /* The generator variables are ordered as [Eqp,Edp,delta,w,Id,Iq] */
235: xgen[idx] = Eqp;
236: xgen[idx+1] = Edp;
237: xgen[idx+2] = delta;
238: xgen[idx+3] = w_s;
240: idx = idx + 4;
242: xgen[idx] = Id;
243: xgen[idx+1] = Iq;
245: idx = idx + 2;
247: /* Exciter */
248: Efd = Eqp + (Xd[i] - Xdp[i])*Id;
249: SE = k1[i]*PetscExpScalar(k2[i]*Efd);
250: VR = KE[i]*Efd + SE;
251: RF = KF[i]*Efd/TF[i];
253: xgen[idx] = Efd;
254: xgen[idx+1] = RF;
255: xgen[idx+2] = VR;
257: Vref[i] = Vm + (VR/KA[i]);
259: idx = idx + 3;
260: }
262: VecRestoreArray(Xgen,&xgen);
263: VecRestoreArray(Xnet,&xnet);
265: /* VecView(Xgen,0); */
266: DMCompositeGather(user->dmpgrid,X,INSERT_VALUES,Xgen,Xnet);
267: DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);
268: return(0);
269: }
273: PetscErrorCode InitialGuess(Vec X,Userctx *user, const PetscScalar PGv[])
274: {
276: Vec Xgen,Xnet;
277: PetscScalar *xgen,*xnet;
278: PetscInt i,idx=0;
279: PetscScalar Vr,Vi,IGr,IGi,Vm,Vm2;
280: PetscScalar Eqp,Edp,delta;
281: PetscScalar Efd,RF,VR; /* Exciter variables */
282: PetscScalar Id,Iq; /* Generator dq axis currents */
283: PetscScalar theta,Vd,Vq,SE;
286: M[0] = 2*H[0]/w_s; M[1] = 2*H[1]/w_s; M[2] = 2*H[2]/w_s;
287: D[0] = 0.1*M[0]; D[1] = 0.1*M[1]; D[2] = 0.1*M[2];
289: DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
291: /* Network subsystem initialization */
292: VecCopy(user->V0,Xnet);
294: /* Generator subsystem initialization */
295: VecGetArray(Xgen,&xgen);
296: VecGetArray(Xnet,&xnet);
298: for (i=0; i < ngen; 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 */
301: Vm = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm;
302: IGr = (Vr*PGv[i] + Vi*QG[i])/Vm2;
303: IGi = (Vi*PGv[i] - Vr*QG[i])/Vm2;
305: delta = atan2(Vi+Xq[i]*IGr,Vr-Xq[i]*IGi); /* Machine angle */
307: theta = PETSC_PI/2.0 - delta;
309: Id = IGr*PetscCosScalar(theta) - IGi*PetscSinScalar(theta); /* d-axis stator current */
310: Iq = IGr*PetscSinScalar(theta) + IGi*PetscCosScalar(theta); /* q-axis stator current */
312: Vd = Vr*PetscCosScalar(theta) - Vi*PetscSinScalar(theta);
313: Vq = Vr*PetscSinScalar(theta) + Vi*PetscCosScalar(theta);
315: Edp = Vd + Rs[i]*Id - Xqp[i]*Iq; /* d-axis transient EMF */
316: Eqp = Vq + Rs[i]*Iq + Xdp[i]*Id; /* q-axis transient EMF */
318: /* The generator variables are ordered as [Eqp,Edp,delta,w,Id,Iq] */
319: xgen[idx] = Eqp;
320: xgen[idx+1] = Edp;
321: xgen[idx+2] = delta;
322: xgen[idx+3] = w_s;
324: idx = idx + 4;
326: xgen[idx] = Id;
327: xgen[idx+1] = Iq;
329: idx = idx + 2;
331: /* Exciter */
332: Efd = Eqp + (Xd[i] - Xdp[i])*Id;
333: SE = k1[i]*PetscExpScalar(k2[i]*Efd);
334: VR = KE[i]*Efd + SE;
335: RF = KF[i]*Efd/TF[i];
337: xgen[idx] = Efd;
338: xgen[idx+1] = RF;
339: xgen[idx+2] = VR;
341: idx = idx + 3;
342: }
344: VecRestoreArray(Xgen,&xgen);
345: VecRestoreArray(Xnet,&xnet);
347: /* VecView(Xgen,0); */
348: DMCompositeGather(user->dmpgrid,X,INSERT_VALUES,Xgen,Xnet);
349: DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);
350: return(0);
351: }
355: PetscErrorCode DICDPFiniteDifference(Vec X,Vec *DICDP, Userctx *user)
356: {
357: Vec Y;
358: PetscScalar PGv[3],eps;
360: PetscInt i,j;
362: eps = 1.e-7;
363: VecDuplicate(X,&Y);
365: for (i=0;i<ngen;i++) {
366: for (j=0;j<3;j++) PGv[j] = PG[j];
367: PGv[i] = PG[i]+eps;
368: InitialGuess(Y,user,PGv);
369: InitialGuess(X,user,PG);
371: VecAXPY(Y,-1.0,X);
372: VecScale(Y,1./eps);
373: VecCopy(Y,DICDP[i]);
374: }
375: VecDestroy(&Y);
376: return(0);
377: }
380: /* Computes F = [-f(x,y);g(x,y)] */
383: PetscErrorCode ResidualFunction(SNES snes,Vec X, Vec F, Userctx *user)
384: {
386: Vec Xgen,Xnet,Fgen,Fnet;
387: PetscScalar *xgen,*xnet,*fgen,*fnet;
388: PetscInt i,idx=0;
389: PetscScalar Vr,Vi,Vm,Vm2;
390: PetscScalar Eqp,Edp,delta,w; /* Generator variables */
391: PetscScalar Efd,RF,VR; /* Exciter variables */
392: PetscScalar Id,Iq; /* Generator dq axis currents */
393: PetscScalar Vd,Vq,SE;
394: PetscScalar IGr,IGi,IDr,IDi;
395: PetscScalar Zdq_inv[4],det;
396: PetscScalar PD,QD,Vm0,*v0;
397: PetscInt k;
400: VecZeroEntries(F);
401: DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
402: DMCompositeGetLocalVectors(user->dmpgrid,&Fgen,&Fnet);
403: DMCompositeScatter(user->dmpgrid,X,Xgen,Xnet);
404: DMCompositeScatter(user->dmpgrid,F,Fgen,Fnet);
406: /* Network current balance residual IG + Y*V + IL = 0. Only YV is added here.
407: The generator current injection, IG, and load current injection, ID are added later
408: */
409: /* Note that the values in Ybus are stored assuming the imaginary current balance
410: equation is ordered first followed by real current balance equation for each bus.
411: Thus imaginary current contribution goes in location 2*i, and
412: real current contribution in 2*i+1
413: */
414: MatMult(user->Ybus,Xnet,Fnet);
416: VecGetArray(Xgen,&xgen);
417: VecGetArray(Xnet,&xnet);
418: VecGetArray(Fgen,&fgen);
419: VecGetArray(Fnet,&fnet);
421: /* Generator subsystem */
422: for (i=0; i < ngen; i++) {
423: Eqp = xgen[idx];
424: Edp = xgen[idx+1];
425: delta = xgen[idx+2];
426: w = xgen[idx+3];
427: Id = xgen[idx+4];
428: Iq = xgen[idx+5];
429: Efd = xgen[idx+6];
430: RF = xgen[idx+7];
431: VR = xgen[idx+8];
433: /* Generator differential equations */
434: fgen[idx] = (Eqp + (Xd[i] - Xdp[i])*Id - Efd)/Td0p[i];
435: fgen[idx+1] = (Edp - (Xq[i] - Xqp[i])*Iq)/Tq0p[i];
436: fgen[idx+2] = -w + w_s;
437: fgen[idx+3] = (-TM[i] + Edp*Id + Eqp*Iq + (Xqp[i] - Xdp[i])*Id*Iq + D[i]*(w - w_s))/M[i];
439: Vr = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */
440: Vi = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */
442: ri2dq(Vr,Vi,delta,&Vd,&Vq);
443: /* Algebraic equations for stator currents */
444: det = Rs[i]*Rs[i] + Xdp[i]*Xqp[i];
446: Zdq_inv[0] = Rs[i]/det;
447: Zdq_inv[1] = Xqp[i]/det;
448: Zdq_inv[2] = -Xdp[i]/det;
449: Zdq_inv[3] = Rs[i]/det;
451: fgen[idx+4] = Zdq_inv[0]*(-Edp + Vd) + Zdq_inv[1]*(-Eqp + Vq) + Id;
452: fgen[idx+5] = Zdq_inv[2]*(-Edp + Vd) + Zdq_inv[3]*(-Eqp + Vq) + Iq;
454: /* Add generator current injection to network */
455: dq2ri(Id,Iq,delta,&IGr,&IGi);
457: fnet[2*gbus[i]] -= IGi;
458: fnet[2*gbus[i]+1] -= IGr;
460: Vm = PetscSqrtScalar(Vd*Vd + Vq*Vq); Vm2 = Vm*Vm;
462: SE = k1[i]*PetscExpScalar(k2[i]*Efd);
464: /* Exciter differential equations */
465: fgen[idx+6] = (KE[i]*Efd + SE - VR)/TE[i];
466: fgen[idx+7] = (RF - KF[i]*Efd/TF[i])/TF[i];
467: fgen[idx+8] = (VR - KA[i]*RF + KA[i]*KF[i]*Efd/TF[i] - KA[i]*(Vref[i] - Vm))/TA[i];
469: idx = idx + 9;
470: }
472: VecGetArray(user->V0,&v0);
473: for (i=0; i < nload; i++) {
474: Vr = xnet[2*lbus[i]]; /* Real part of load bus voltage */
475: Vi = xnet[2*lbus[i]+1]; /* Imaginary part of the load bus voltage */
476: Vm = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm;
477: Vm0 = PetscSqrtScalar(v0[2*lbus[i]]*v0[2*lbus[i]] + v0[2*lbus[i]+1]*v0[2*lbus[i]+1]);
478: PD = QD = 0.0;
479: for (k=0; k < ld_nsegsp[i]; k++) PD += ld_alphap[k]*PD0[i]*PetscPowScalar((Vm/Vm0),ld_betap[k]);
480: for (k=0; k < ld_nsegsq[i]; k++) QD += ld_alphaq[k]*QD0[i]*PetscPowScalar((Vm/Vm0),ld_betaq[k]);
482: /* Load currents */
483: IDr = (PD*Vr + QD*Vi)/Vm2;
484: IDi = (-QD*Vr + PD*Vi)/Vm2;
486: fnet[2*lbus[i]] += IDi;
487: fnet[2*lbus[i]+1] += IDr;
488: }
489: VecRestoreArray(user->V0,&v0);
491: VecRestoreArray(Xgen,&xgen);
492: VecRestoreArray(Xnet,&xnet);
493: VecRestoreArray(Fgen,&fgen);
494: VecRestoreArray(Fnet,&fnet);
496: DMCompositeGather(user->dmpgrid,F,INSERT_VALUES,Fgen,Fnet);
497: DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);
498: DMCompositeRestoreLocalVectors(user->dmpgrid,&Fgen,&Fnet);
499: return(0);
500: }
502: /* \dot{x} - f(x,y)
503: g(x,y) = 0
504: */
507: PetscErrorCode IFunction(TS ts,PetscReal t, Vec X, Vec Xdot, Vec F, Userctx *user)
508: {
509: PetscErrorCode ierr;
510: SNES snes;
511: PetscScalar *f;
512: const PetscScalar *xdot;
513: PetscInt i;
516: user->t = t;
518: TSGetSNES(ts,&snes);
519: ResidualFunction(snes,X,F,user);
520: VecGetArray(F,&f);
521: VecGetArrayRead(Xdot,&xdot);
522: for (i=0;i < ngen;i++) {
523: f[9*i] += xdot[9*i];
524: f[9*i+1] += xdot[9*i+1];
525: f[9*i+2] += xdot[9*i+2];
526: f[9*i+3] += xdot[9*i+3];
527: f[9*i+6] += xdot[9*i+6];
528: f[9*i+7] += xdot[9*i+7];
529: f[9*i+8] += xdot[9*i+8];
530: }
531: VecRestoreArray(F,&f);
532: VecRestoreArrayRead(Xdot,&xdot);
533: return(0);
534: }
536: /* This function is used for solving the algebraic system only during fault on and
537: off times. It computes the entire F and then zeros out the part corresponding to
538: differential equations
539: F = [0;g(y)];
540: */
543: PetscErrorCode AlgFunction(SNES snes, Vec X, Vec F, void *ctx)
544: {
546: Userctx *user=(Userctx*)ctx;
547: PetscScalar *f;
548: PetscInt i;
551: ResidualFunction(snes,X,F,user);
552: VecGetArray(F,&f);
553: for (i=0; i < ngen; i++) {
554: f[9*i] = 0;
555: f[9*i+1] = 0;
556: f[9*i+2] = 0;
557: f[9*i+3] = 0;
558: f[9*i+6] = 0;
559: f[9*i+7] = 0;
560: f[9*i+8] = 0;
561: }
562: VecRestoreArray(F,&f);
563: return(0);
564: }
568: PetscErrorCode PreallocateJacobian(Mat J, Userctx *user)
569: {
571: PetscInt *d_nnz;
572: PetscInt i,idx=0,start=0;
573: PetscInt ncols;
576: PetscMalloc1(user->neqs_pgrid,&d_nnz);
577: for (i=0; i<user->neqs_pgrid; i++) d_nnz[i] = 0;
578: /* Generator subsystem */
579: for (i=0; i < ngen; i++) {
581: d_nnz[idx] += 3;
582: d_nnz[idx+1] += 2;
583: d_nnz[idx+2] += 2;
584: d_nnz[idx+3] += 5;
585: d_nnz[idx+4] += 6;
586: d_nnz[idx+5] += 6;
588: d_nnz[user->neqs_gen+2*gbus[i]] += 3;
589: d_nnz[user->neqs_gen+2*gbus[i]+1] += 3;
591: d_nnz[idx+6] += 2;
592: d_nnz[idx+7] += 2;
593: d_nnz[idx+8] += 5;
595: idx = idx + 9;
596: }
598: start = user->neqs_gen;
599: for (i=0; i < nbus; i++) {
600: MatGetRow(user->Ybus,2*i,&ncols,NULL,NULL);
601: d_nnz[start+2*i] += ncols;
602: d_nnz[start+2*i+1] += ncols;
603: MatRestoreRow(user->Ybus,2*i,&ncols,NULL,NULL);
604: }
606: MatSeqAIJSetPreallocation(J,0,d_nnz);
607: PetscFree(d_nnz);
608: return(0);
609: }
611: /*
612: J = [-df_dx, -df_dy
613: dg_dx, dg_dy]
614: */
617: PetscErrorCode ResidualJacobian(SNES snes,Vec X,Mat J,Mat B,void *ctx)
618: {
619: PetscErrorCode ierr;
620: Userctx *user=(Userctx*)ctx;
621: Vec Xgen,Xnet;
622: PetscScalar *xgen,*xnet;
623: PetscInt i,idx=0;
624: PetscScalar Vr,Vi,Vm,Vm2;
625: PetscScalar Eqp,Edp,delta; /* Generator variables */
626: PetscScalar Efd; /* Exciter variables */
627: PetscScalar Id,Iq; /* Generator dq axis currents */
628: PetscScalar Vd,Vq;
629: PetscScalar val[10];
630: PetscInt row[2],col[10];
631: PetscInt net_start=user->neqs_gen;
632: PetscInt ncols;
633: const PetscInt *cols;
634: const PetscScalar *yvals;
635: PetscInt k;
636: PetscScalar Zdq_inv[4],det;
637: PetscScalar dVd_dVr,dVd_dVi,dVq_dVr,dVq_dVi,dVd_ddelta,dVq_ddelta;
638: PetscScalar dIGr_ddelta,dIGi_ddelta,dIGr_dId,dIGr_dIq,dIGi_dId,dIGi_dIq;
639: PetscScalar dSE_dEfd;
640: PetscScalar dVm_dVd,dVm_dVq,dVm_dVr,dVm_dVi;
641: PetscScalar PD,QD,Vm0,*v0,Vm4;
642: PetscScalar dPD_dVr,dPD_dVi,dQD_dVr,dQD_dVi;
643: PetscScalar dIDr_dVr,dIDr_dVi,dIDi_dVr,dIDi_dVi;
646: MatZeroEntries(B);
647: DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
648: DMCompositeScatter(user->dmpgrid,X,Xgen,Xnet);
650: VecGetArray(Xgen,&xgen);
651: VecGetArray(Xnet,&xnet);
653: /* Generator subsystem */
654: for (i=0; i < ngen; i++) {
655: Eqp = xgen[idx];
656: Edp = xgen[idx+1];
657: delta = xgen[idx+2];
658: Id = xgen[idx+4];
659: Iq = xgen[idx+5];
660: Efd = xgen[idx+6];
662: /* fgen[idx] = (Eqp + (Xd[i] - Xdp[i])*Id - Efd)/Td0p[i]; */
663: row[0] = idx;
664: col[0] = idx; col[1] = idx+4; col[2] = idx+6;
665: val[0] = 1/ Td0p[i]; val[1] = (Xd[i] - Xdp[i])/ Td0p[i]; val[2] = -1/Td0p[i];
667: MatSetValues(J,1,row,3,col,val,INSERT_VALUES);
669: /* fgen[idx+1] = (Edp - (Xq[i] - Xqp[i])*Iq)/Tq0p[i]; */
670: row[0] = idx + 1;
671: col[0] = idx + 1; col[1] = idx+5;
672: val[0] = 1/Tq0p[i]; val[1] = -(Xq[i] - Xqp[i])/Tq0p[i];
673: MatSetValues(J,1,row,2,col,val,INSERT_VALUES);
675: /* fgen[idx+2] = - w + w_s; */
676: row[0] = idx + 2;
677: col[0] = idx + 2; col[1] = idx + 3;
678: val[0] = 0; val[1] = -1;
679: MatSetValues(J,1,row,2,col,val,INSERT_VALUES);
681: /* fgen[idx+3] = (-TM[i] + Edp*Id + Eqp*Iq + (Xqp[i] - Xdp[i])*Id*Iq + D[i]*(w - w_s))/M[i]; */
682: row[0] = idx + 3;
683: col[0] = idx; col[1] = idx + 1; col[2] = idx + 3; col[3] = idx + 4; col[4] = idx + 5;
684: 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];
685: MatSetValues(J,1,row,5,col,val,INSERT_VALUES);
687: Vr = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */
688: Vi = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */
689: ri2dq(Vr,Vi,delta,&Vd,&Vq);
691: det = Rs[i]*Rs[i] + Xdp[i]*Xqp[i];
693: Zdq_inv[0] = Rs[i]/det;
694: Zdq_inv[1] = Xqp[i]/det;
695: Zdq_inv[2] = -Xdp[i]/det;
696: Zdq_inv[3] = Rs[i]/det;
698: dVd_dVr = PetscSinScalar(delta); dVd_dVi = -PetscCosScalar(delta);
699: dVq_dVr = PetscCosScalar(delta); dVq_dVi = PetscSinScalar(delta);
700: dVd_ddelta = Vr*PetscCosScalar(delta) + Vi*PetscSinScalar(delta);
701: dVq_ddelta = -Vr*PetscSinScalar(delta) + Vi*PetscCosScalar(delta);
703: /* fgen[idx+4] = Zdq_inv[0]*(-Edp + Vd) + Zdq_inv[1]*(-Eqp + Vq) + Id; */
704: row[0] = idx+4;
705: col[0] = idx; col[1] = idx+1; col[2] = idx + 2;
706: val[0] = -Zdq_inv[1]; val[1] = -Zdq_inv[0]; val[2] = Zdq_inv[0]*dVd_ddelta + Zdq_inv[1]*dVq_ddelta;
707: col[3] = idx + 4; col[4] = net_start+2*gbus[i]; col[5] = net_start + 2*gbus[i]+1;
708: 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;
709: MatSetValues(J,1,row,6,col,val,INSERT_VALUES);
711: /* fgen[idx+5] = Zdq_inv[2]*(-Edp + Vd) + Zdq_inv[3]*(-Eqp + Vq) + Iq; */
712: row[0] = idx+5;
713: col[0] = idx; col[1] = idx+1; col[2] = idx + 2;
714: val[0] = -Zdq_inv[3]; val[1] = -Zdq_inv[2]; val[2] = Zdq_inv[2]*dVd_ddelta + Zdq_inv[3]*dVq_ddelta;
715: col[3] = idx + 5; col[4] = net_start+2*gbus[i]; col[5] = net_start + 2*gbus[i]+1;
716: 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;
717: MatSetValues(J,1,row,6,col,val,INSERT_VALUES);
719: dIGr_ddelta = Id*PetscCosScalar(delta) - Iq*PetscSinScalar(delta);
720: dIGi_ddelta = Id*PetscSinScalar(delta) + Iq*PetscCosScalar(delta);
721: dIGr_dId = PetscSinScalar(delta); dIGr_dIq = PetscCosScalar(delta);
722: dIGi_dId = -PetscCosScalar(delta); dIGi_dIq = PetscSinScalar(delta);
724: /* fnet[2*gbus[i]] -= IGi; */
725: row[0] = net_start + 2*gbus[i];
726: col[0] = idx+2; col[1] = idx + 4; col[2] = idx + 5;
727: val[0] = -dIGi_ddelta; val[1] = -dIGi_dId; val[2] = -dIGi_dIq;
728: MatSetValues(J,1,row,3,col,val,INSERT_VALUES);
730: /* fnet[2*gbus[i]+1] -= IGr; */
731: row[0] = net_start + 2*gbus[i]+1;
732: col[0] = idx+2; col[1] = idx + 4; col[2] = idx + 5;
733: val[0] = -dIGr_ddelta; val[1] = -dIGr_dId; val[2] = -dIGr_dIq;
734: MatSetValues(J,1,row,3,col,val,INSERT_VALUES);
736: Vm = PetscSqrtScalar(Vd*Vd + Vq*Vq); Vm2 = Vm*Vm;
738: /* fgen[idx+6] = (KE[i]*Efd + SE - VR)/TE[i]; */
739: /* SE = k1[i]*PetscExpScalar(k2[i]*Efd); */
740: dSE_dEfd = k1[i]*k2[i]*PetscExpScalar(k2[i]*Efd);
742: row[0] = idx + 6;
743: col[0] = idx + 6; col[1] = idx + 8;
744: val[0] = (KE[i] + dSE_dEfd)/TE[i]; val[1] = -1/TE[i];
745: MatSetValues(J,1,row,2,col,val,INSERT_VALUES);
747: /* Exciter differential equations */
749: /* fgen[idx+7] = (RF - KF[i]*Efd/TF[i])/TF[i]; */
750: row[0] = idx + 7;
751: col[0] = idx + 6; col[1] = idx + 7;
752: val[0] = (-KF[i]/TF[i])/TF[i]; val[1] = 1/TF[i];
753: MatSetValues(J,1,row,2,col,val,INSERT_VALUES);
755: /* fgen[idx+8] = (VR - KA[i]*RF + KA[i]*KF[i]*Efd/TF[i] - KA[i]*(Vref[i] - Vm))/TA[i]; */
756: /* Vm = (Vd^2 + Vq^2)^0.5; */
757: dVm_dVd = Vd/Vm; dVm_dVq = Vq/Vm;
758: dVm_dVr = dVm_dVd*dVd_dVr + dVm_dVq*dVq_dVr;
759: dVm_dVi = dVm_dVd*dVd_dVi + dVm_dVq*dVq_dVi;
760: row[0] = idx + 8;
761: col[0] = idx + 6; col[1] = idx + 7; col[2] = idx + 8;
762: val[0] = (KA[i]*KF[i]/TF[i])/TA[i]; val[1] = -KA[i]/TA[i]; val[2] = 1/TA[i];
763: col[3] = net_start + 2*gbus[i]; col[4] = net_start + 2*gbus[i]+1;
764: val[3] = KA[i]*dVm_dVr/TA[i]; val[4] = KA[i]*dVm_dVi/TA[i];
765: MatSetValues(J,1,row,5,col,val,INSERT_VALUES);
766: idx = idx + 9;
767: }
770: for (i=0; i<nbus; i++) {
771: MatGetRow(user->Ybus,2*i,&ncols,&cols,&yvals);
772: row[0] = net_start + 2*i;
773: for (k=0; k<ncols; k++) {
774: col[k] = net_start + cols[k];
775: val[k] = yvals[k];
776: }
777: MatSetValues(J,1,row,ncols,col,val,INSERT_VALUES);
778: MatRestoreRow(user->Ybus,2*i,&ncols,&cols,&yvals);
780: MatGetRow(user->Ybus,2*i+1,&ncols,&cols,&yvals);
781: row[0] = net_start + 2*i+1;
782: for (k=0; k<ncols; k++) {
783: col[k] = net_start + cols[k];
784: val[k] = yvals[k];
785: }
786: MatSetValues(J,1,row,ncols,col,val,INSERT_VALUES);
787: MatRestoreRow(user->Ybus,2*i+1,&ncols,&cols,&yvals);
788: }
790: MatAssemblyBegin(J,MAT_FLUSH_ASSEMBLY);
791: MatAssemblyEnd(J,MAT_FLUSH_ASSEMBLY);
794: VecGetArray(user->V0,&v0);
795: for (i=0; i < nload; i++) {
796: Vr = xnet[2*lbus[i]]; /* Real part of load bus voltage */
797: Vi = xnet[2*lbus[i]+1]; /* Imaginary part of the load bus voltage */
798: Vm = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm; Vm4 = Vm2*Vm2;
799: Vm0 = PetscSqrtScalar(v0[2*lbus[i]]*v0[2*lbus[i]] + v0[2*lbus[i]+1]*v0[2*lbus[i]+1]);
800: PD = QD = 0.0;
801: dPD_dVr = dPD_dVi = dQD_dVr = dQD_dVi = 0.0;
802: for (k=0; k < ld_nsegsp[i]; k++) {
803: PD += ld_alphap[k]*PD0[i]*PetscPowScalar((Vm/Vm0),ld_betap[k]);
804: dPD_dVr += ld_alphap[k]*ld_betap[k]*PD0[i]*PetscPowScalar((1/Vm0),ld_betap[k])*Vr*PetscPowScalar(Vm,(ld_betap[k]-2));
805: dPD_dVi += ld_alphap[k]*ld_betap[k]*PD0[i]*PetscPowScalar((1/Vm0),ld_betap[k])*Vi*PetscPowScalar(Vm,(ld_betap[k]-2));
806: }
807: for (k=0; k < ld_nsegsq[i]; k++) {
808: QD += ld_alphaq[k]*QD0[i]*PetscPowScalar((Vm/Vm0),ld_betaq[k]);
809: dQD_dVr += ld_alphaq[k]*ld_betaq[k]*QD0[i]*PetscPowScalar((1/Vm0),ld_betaq[k])*Vr*PetscPowScalar(Vm,(ld_betaq[k]-2));
810: dQD_dVi += ld_alphaq[k]*ld_betaq[k]*QD0[i]*PetscPowScalar((1/Vm0),ld_betaq[k])*Vi*PetscPowScalar(Vm,(ld_betaq[k]-2));
811: }
813: /* IDr = (PD*Vr + QD*Vi)/Vm2; */
814: /* IDi = (-QD*Vr + PD*Vi)/Vm2; */
816: dIDr_dVr = (dPD_dVr*Vr + dQD_dVr*Vi + PD)/Vm2 - ((PD*Vr + QD*Vi)*2*Vr)/Vm4;
817: dIDr_dVi = (dPD_dVi*Vr + dQD_dVi*Vi + QD)/Vm2 - ((PD*Vr + QD*Vi)*2*Vi)/Vm4;
819: dIDi_dVr = (-dQD_dVr*Vr + dPD_dVr*Vi - QD)/Vm2 - ((-QD*Vr + PD*Vi)*2*Vr)/Vm4;
820: dIDi_dVi = (-dQD_dVi*Vr + dPD_dVi*Vi + PD)/Vm2 - ((-QD*Vr + PD*Vi)*2*Vi)/Vm4;
823: /* fnet[2*lbus[i]] += IDi; */
824: row[0] = net_start + 2*lbus[i];
825: col[0] = net_start + 2*lbus[i]; col[1] = net_start + 2*lbus[i]+1;
826: val[0] = dIDi_dVr; val[1] = dIDi_dVi;
827: MatSetValues(J,1,row,2,col,val,ADD_VALUES);
828: /* fnet[2*lbus[i]+1] += IDr; */
829: row[0] = net_start + 2*lbus[i]+1;
830: col[0] = net_start + 2*lbus[i]; col[1] = net_start + 2*lbus[i]+1;
831: val[0] = dIDr_dVr; val[1] = dIDr_dVi;
832: MatSetValues(J,1,row,2,col,val,ADD_VALUES);
833: }
834: VecRestoreArray(user->V0,&v0);
836: VecRestoreArray(Xgen,&xgen);
837: VecRestoreArray(Xnet,&xnet);
839: DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);
841: MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
842: MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
843: return(0);
844: }
846: /*
847: J = [I, 0
848: dg_dx, dg_dy]
849: */
852: PetscErrorCode AlgJacobian(SNES snes,Vec X,Mat A,Mat B,void *ctx)
853: {
855: Userctx *user=(Userctx*)ctx;
858: ResidualJacobian(snes,X,A,B,ctx);
859: MatSetOption(A,MAT_KEEP_NONZERO_PATTERN,PETSC_TRUE);
860: MatZeroRowsIS(A,user->is_diff,1.0,NULL,NULL);
861: return(0);
862: }
864: /*
865: J = [a*I-df_dx, -df_dy
866: dg_dx, dg_dy]
867: */
871: PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,Userctx *user)
872: {
874: SNES snes;
875: PetscScalar atmp = (PetscScalar) a;
876: PetscInt i,row;
879: user->t = t;
881: TSGetSNES(ts,&snes);
882: ResidualJacobian(snes,X,A,B,user);
883: for (i=0;i < ngen;i++) {
884: row = 9*i;
885: MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
886: row = 9*i+1;
887: MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
888: row = 9*i+2;
889: MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
890: row = 9*i+3;
891: MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
892: row = 9*i+6;
893: MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
894: row = 9*i+7;
895: MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
896: row = 9*i+8;
897: MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
898: }
899: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
900: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
901: return(0);
902: }
904: /* Matrix JacobianP is constant so that it only needs to be evaluated once */
907: static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A, void *ctx0)
908: {
910: PetscScalar a;
911: PetscInt row,col;
912: Userctx *ctx=(Userctx*)ctx0;
916: if (ctx->jacp_flg) {
917: MatZeroEntries(A);
919: for (col=0;col<3;col++) {
920: a = 1.0/M[col];
921: row = 9*col+3;
922: MatSetValues(A,1,&row,1,&col,&a,INSERT_VALUES);
923: }
925: ctx->jacp_flg = PETSC_FALSE;
926:
927: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
928: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
929: }
930: return(0);
931: }
935: static PetscErrorCode CostIntegrand(TS ts,PetscReal t,Vec U,Vec R,Userctx *user)
936: {
938: PetscScalar *u,*r;
939: PetscInt idx;
940: Vec Xgen,Xnet;
941: PetscScalar *xgen;
942: PetscInt i;
945: DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
946: DMCompositeScatter(user->dmpgrid,U,Xgen,Xnet);
948: VecGetArray(Xgen,&xgen);
950: VecGetArray(U,&u);
951: VecGetArray(R,&r);
952: r[0] = 0.;
953: idx = 0;
954: for (i=0;i<ngen;i++) {
955: 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);
956: idx += 9;
957: }
958: VecRestoreArray(R,&r);
959: DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);
960: return(0);
961: }
965: static PetscErrorCode DRDYFunction(TS ts,PetscReal t,Vec U,Vec *drdy,Userctx *user)
966: {
968: Vec Xgen,Xnet,Dgen,Dnet;
969: PetscScalar *xgen,*dgen;
970: PetscInt i;
971: PetscInt idx;
974: DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
975: DMCompositeGetLocalVectors(user->dmpgrid,&Dgen,&Dnet);
976: DMCompositeScatter(user->dmpgrid,U,Xgen,Xnet);
977: DMCompositeScatter(user->dmpgrid,drdy[0],Dgen,Dnet);
979: VecGetArray(Xgen,&xgen);
980: VecGetArray(Dgen,&dgen);
982: idx = 0;
983: for (i=0;i<ngen;i++) {
984: dgen[idx+3] = 0.;
985: if (xgen[idx+3]/(2.*PETSC_PI) > user->freq_u) dgen[idx+3] = user->pow*PetscPowScalarInt(xgen[idx+3]/(2.*PETSC_PI)-user->freq_u,user->pow-1)/(2.*PETSC_PI);
986: if (xgen[idx+3]/(2.*PETSC_PI) < user->freq_l) dgen[idx+3] = user->pow*PetscPowScalarInt(user->freq_l-xgen[idx+3]/(2.*PETSC_PI),user->pow-1)/(-2.*PETSC_PI);
987: idx += 9;
988: }
990: VecRestoreArray(Dgen,&dgen);
991: DMCompositeGather(user->dmpgrid,drdy[0],INSERT_VALUES,Dgen,Dnet);
992: DMCompositeRestoreLocalVectors(user->dmpgrid,&Dgen,&Dnet);
993: DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);
994: return(0);
995: }
999: static PetscErrorCode DRDPFunction(TS ts,PetscReal t,Vec U,Vec *drdp,Userctx *user)
1000: {
1002: return(0);
1003: }
1007: PetscErrorCode ComputeSensiP(Vec lambda,Vec mu,Vec *DICDP,Userctx *user)
1008: {
1010: PetscScalar *x,*y,sensip;
1011: PetscInt i;
1014: VecGetArray(lambda,&x);
1015: VecGetArray(mu,&y);
1017: for (i=0;i<3;i++) {
1018: VecDot(lambda,DICDP[i],&sensip);
1019: sensip = sensip+y[i];
1020: /* PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt %D th parameter: %g \n",i,(double)sensip); */
1021: y[i] = sensip;
1022: }
1023: VecRestoreArray(mu,&y);
1024: return(0);
1025: }
1029: int main(int argc,char **argv)
1030: {
1031: Userctx user;
1032: Vec p;
1033: PetscScalar *x_ptr;
1034: PetscErrorCode ierr;
1035: PetscMPIInt size;
1036: PetscInt i;
1037: PetscViewer Xview,Ybusview;
1038: PetscInt *idx2;
1039: Tao tao;
1040: TaoConvergedReason reason;
1041: KSP ksp;
1042: PC pc;
1043: Vec lowerb,upperb;
1045: PetscInitialize(&argc,&argv,"petscoptions",help);
1046: MPI_Comm_size(PETSC_COMM_WORLD,&size);
1047: if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs");
1049: user.jacp_flg = PETSC_TRUE;
1050: user.neqs_gen = 9*ngen; /* # eqs. for generator subsystem */
1051: user.neqs_net = 2*nbus; /* # eqs. for network subsystem */
1052: user.neqs_pgrid = user.neqs_gen + user.neqs_net;
1054: /* Create indices for differential and algebraic equations */
1055: PetscMalloc1(7*ngen,&idx2);
1056: for (i=0; i<ngen; i++) {
1057: 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;
1058: idx2[7*i+4] = 9*i+6; idx2[7*i+5] = 9*i+7; idx2[7*i+6] = 9*i+8;
1059: }
1060: ISCreateGeneral(PETSC_COMM_WORLD,7*ngen,idx2,PETSC_COPY_VALUES,&user.is_diff);
1061: ISComplement(user.is_diff,0,user.neqs_pgrid,&user.is_alg);
1062: PetscFree(idx2);
1064: /* Set run time options */
1065: PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Transient stability fault options","");
1066: {
1067: user.tfaulton = 1.0;
1068: user.tfaultoff = 1.2;
1069: user.Rfault = 0.0001;
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: user.freq_u = 61.0;
1079: user.freq_l = 59.0;
1080: user.pow = 2;
1081: PetscOptionsReal("-frequ","","",user.freq_u,&user.freq_u,NULL);
1082: PetscOptionsReal("-freql","","",user.freq_l,&user.freq_l,NULL);
1083: PetscOptionsInt("-pow","","",user.pow,&user.pow,NULL);
1085: }
1086: PetscOptionsEnd();
1088: /* Create DMs for generator and network subsystems */
1089: DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,user.neqs_gen,1,1,NULL,&user.dmgen);
1090: DMSetOptionsPrefix(user.dmgen,"dmgen_");
1091: DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,user.neqs_net,1,1,NULL,&user.dmnet);
1092: DMSetOptionsPrefix(user.dmnet,"dmnet_");
1093: /* Create a composite DM packer and add the two DMs */
1094: DMCompositeCreate(PETSC_COMM_WORLD,&user.dmpgrid);
1095: DMSetOptionsPrefix(user.dmpgrid,"pgrid_");
1096: DMCompositeAddDM(user.dmpgrid,user.dmgen);
1097: DMCompositeAddDM(user.dmpgrid,user.dmnet);
1099: /* Read initial voltage vector and Ybus */
1100: PetscViewerBinaryOpen(PETSC_COMM_WORLD,"X.bin",FILE_MODE_READ,&Xview);
1101: PetscViewerBinaryOpen(PETSC_COMM_WORLD,"Ybus.bin",FILE_MODE_READ,&Ybusview);
1103: VecCreate(PETSC_COMM_WORLD,&user.V0);
1104: VecSetSizes(user.V0,PETSC_DECIDE,user.neqs_net);
1105: VecLoad(user.V0,Xview);
1107: MatCreate(PETSC_COMM_WORLD,&user.Ybus);
1108: MatSetSizes(user.Ybus,PETSC_DECIDE,PETSC_DECIDE,user.neqs_net,user.neqs_net);
1109: MatSetType(user.Ybus,MATBAIJ);
1110: /* MatSetBlockSize(ctx->Ybus,2); */
1111: MatLoad(user.Ybus,Ybusview);
1113: PetscViewerDestroy(&Xview);
1114: PetscViewerDestroy(&Ybusview);
1116: /* Allocate space for Jacobians */
1117: MatCreate(PETSC_COMM_WORLD,&user.J);
1118: MatSetSizes(user.J,PETSC_DECIDE,PETSC_DECIDE,user.neqs_pgrid,user.neqs_pgrid);
1119: MatSetFromOptions(user.J);
1120: PreallocateJacobian(user.J,&user);
1122: MatCreate(PETSC_COMM_WORLD,&user.Jacp);
1123: MatSetSizes(user.Jacp,PETSC_DECIDE,PETSC_DECIDE,user.neqs_pgrid,3);
1124: MatSetFromOptions(user.Jacp);
1125: MatSetUp(user.Jacp);
1126: MatZeroEntries(user.Jacp); /* initialize to zeros */
1128: /* Create TAO solver and set desired solution method */
1129: TaoCreate(PETSC_COMM_WORLD,&tao);
1130: TaoSetType(tao,TAOBLMVM);
1131: /*
1132: Optimization starts
1133: */
1134: /* Set initial solution guess */
1135: VecCreateSeq(PETSC_COMM_WORLD,3,&p);
1136: VecGetArray(p,&x_ptr);
1137: x_ptr[0] = PG[0]; x_ptr[1] = PG[1]; x_ptr[2] = PG[2];
1138: VecRestoreArray(p,&x_ptr);
1140: TaoSetInitialVector(tao,p);
1141: /* Set routine for function and gradient evaluation */
1142: TaoSetObjectiveAndGradientRoutine(tao,FormFunctionGradient,&user);
1143:
1144: /* Set bounds for the optimization */
1145: VecDuplicate(p,&lowerb);
1146: VecDuplicate(p,&upperb);
1147: VecGetArray(lowerb,&x_ptr);
1148: x_ptr[0] = 0.5; x_ptr[1] = 0.5; x_ptr[2] = 0.5;
1149: VecRestoreArray(lowerb,&x_ptr);
1150: VecGetArray(upperb,&x_ptr);
1151: x_ptr[0] = 2.0; x_ptr[1] = 2.0; x_ptr[2] = 2.0;
1152: VecRestoreArray(upperb,&x_ptr);
1153: TaoSetVariableBounds(tao,lowerb,upperb);
1155: /* Check for any TAO command line options */
1156: TaoSetFromOptions(tao);
1157: TaoGetKSP(tao,&ksp);
1158: if (ksp) {
1159: KSPGetPC(ksp,&pc);
1160: PCSetType(pc,PCNONE);
1161: }
1163: /* TaoSetTolerances(tao,1e-15,1e-15,1e-15,1e-15,1e-15); */
1164: /* SOLVE THE APPLICATION */
1165: TaoSolve(tao);
1166: /* Get information on termination */
1167: TaoGetConvergedReason(tao,&reason);
1168: if (reason <= 0){
1169: ierr=PetscPrintf(MPI_COMM_WORLD, "Try another method! \n");
1170: }
1172: VecView(p,PETSC_VIEWER_STDOUT_WORLD);
1173: /* Free TAO data structures */
1174: TaoDestroy(&tao);
1176: DMDestroy(&user.dmgen);
1177: DMDestroy(&user.dmnet);
1178: DMDestroy(&user.dmpgrid);
1179: ISDestroy(&user.is_diff);
1180: ISDestroy(&user.is_alg);
1181:
1182: MatDestroy(&user.J);
1183: MatDestroy(&user.Jacp);
1184: MatDestroy(&user.Ybus);
1185: /* MatDestroy(&user.Sol); */
1186: VecDestroy(&user.V0);
1187: VecDestroy(&p);
1188: VecDestroy(&lowerb);
1189: VecDestroy(&upperb);
1190: PetscFinalize();
1191: return(0);
1192: }
1194: /* ------------------------------------------------------------------ */
1197: /*
1198: FormFunction - Evaluates the function and corresponding gradient.
1200: Input Parameters:
1201: tao - the Tao context
1202: X - the input vector
1203: ptr - optional user-defined context, as set by TaoSetObjectiveAndGradientRoutine()
1205: Output Parameters:
1206: f - the newly evaluated function
1207: G - the newly evaluated gradient
1208: */
1209: PetscErrorCode FormFunctionGradient(Tao tao,Vec P,PetscReal *f,Vec G,void *ctx0)
1210: {
1211: TS ts;
1212: SNES snes_alg;
1214: Userctx *ctx = (Userctx*)ctx0;
1215: Vec X;
1216: PetscInt i;
1217: /* sensitivity context */
1218: PetscScalar *x_ptr;
1219: Vec lambda[1],q;
1220: Vec mu[1];
1221: PetscInt steps1,steps2,steps3;
1222: Vec DICDP[3];
1223: Vec F_alg;
1224: PetscInt row_loc,col_loc;
1225: PetscScalar val;
1226: Vec Xdot;
1228: VecGetArray(P,&x_ptr);
1229: PG[0] = x_ptr[0];
1230: PG[1] = x_ptr[1];
1231: PG[2] = x_ptr[2];
1232: VecRestoreArray(P,&x_ptr);
1234: ctx->stepnum = 0;
1236: DMCreateGlobalVector(ctx->dmpgrid,&X);
1238: /* Create matrix to save solutions at each time step */
1239: /* MatCreateSeqDense(PETSC_COMM_SELF,ctx->neqs_pgrid+1,1002,NULL,&ctx->Sol); */
1240: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1241: Create timestepping solver context
1242: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1243: TSCreate(PETSC_COMM_WORLD,&ts);
1244: TSSetProblemType(ts,TS_NONLINEAR);
1245: TSSetType(ts,TSCN);
1246: TSSetIFunction(ts,NULL,(TSIFunction) IFunction,ctx);
1247: TSSetIJacobian(ts,ctx->J,ctx->J,(TSIJacobian)IJacobian,ctx);
1248: TSSetApplicationContext(ts,ctx);
1250: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1251: Set initial conditions
1252: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1253: SetInitialGuess(X,ctx);
1255: /* Approximate DICDP with finite difference, we want to zero out network variables */
1256: for (i=0;i<3;i++) {
1257: VecDuplicate(X,&DICDP[i]);
1258: }
1259: DICDPFiniteDifference(X,DICDP,ctx);
1261: VecDuplicate(X,&F_alg);
1262: SNESCreate(PETSC_COMM_WORLD,&snes_alg);
1263: SNESSetFunction(snes_alg,F_alg,AlgFunction,ctx);
1264: MatZeroEntries(ctx->J);
1265: SNESSetJacobian(snes_alg,ctx->J,ctx->J,AlgJacobian,ctx);
1266: SNESSetOptionsPrefix(snes_alg,"alg_");
1267: SNESSetFromOptions(snes_alg);
1268: ctx->alg_flg = PETSC_TRUE;
1269: /* Solve the algebraic equations */
1270: SNESSolve(snes_alg,NULL,X);
1272: /* Just to set up the Jacobian structure */
1273: VecDuplicate(X,&Xdot);
1274: IJacobian(ts,0.0,X,Xdot,0.0,ctx->J,ctx->J,ctx);
1275: VecDestroy(&Xdot);
1277: ctx->stepnum++;
1279: /*
1280: Save trajectory of solution so that TSAdjointSolve() may be used
1281: */
1282: TSSetSaveTrajectory(ts);
1284: TSSetDuration(ts,1000,ctx->tfaulton);
1285: TSSetInitialTimeStep(ts,0.0,0.01);
1286: TSSetFromOptions(ts);
1287: /* TSSetPostStep(ts,SaveSolution); */
1289: ctx->alg_flg = PETSC_FALSE;
1290: /* Prefault period */
1291: TSSolve(ts,X);
1292: TSGetTimeStepNumber(ts,&steps1);
1294: /* Create the nonlinear solver for solving the algebraic system */
1295: /* Note that although the algebraic system needs to be solved only for
1296: Idq and V, we reuse the entire system including xgen. The xgen
1297: variables are held constant by setting their residuals to 0 and
1298: putting a 1 on the Jacobian diagonal for xgen rows
1299: */
1300: MatZeroEntries(ctx->J);
1302: /* Apply disturbance - resistive fault at ctx->faultbus */
1303: /* This is done by adding shunt conductance to the diagonal location
1304: in the Ybus matrix */
1305: row_loc = 2*ctx->faultbus; col_loc = 2*ctx->faultbus+1; /* Location for G */
1306: val = 1/ctx->Rfault;
1307: MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);
1308: row_loc = 2*ctx->faultbus+1; col_loc = 2*ctx->faultbus; /* Location for G */
1309: val = 1/ctx->Rfault;
1310: MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);
1312: MatAssemblyBegin(ctx->Ybus,MAT_FINAL_ASSEMBLY);
1313: MatAssemblyEnd(ctx->Ybus,MAT_FINAL_ASSEMBLY);
1315: ctx->alg_flg = PETSC_TRUE;
1316: /* Solve the algebraic equations */
1317: SNESSolve(snes_alg,NULL,X);
1319: ctx->stepnum++;
1321: /* Disturbance period */
1322: TSSetDuration(ts,1000,ctx->tfaultoff);
1323: TSSetInitialTimeStep(ts,ctx->tfaulton,.01);
1325: ctx->alg_flg = PETSC_FALSE;
1327: TSSolve(ts,X);
1328: TSGetTimeStepNumber(ts,&steps2);
1330: /* Remove the fault */
1331: row_loc = 2*ctx->faultbus; col_loc = 2*ctx->faultbus+1;
1332: val = -1/ctx->Rfault;
1333: MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);
1334: row_loc = 2*ctx->faultbus+1; col_loc = 2*ctx->faultbus;
1335: val = -1/ctx->Rfault;
1336: MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);
1338: MatAssemblyBegin(ctx->Ybus,MAT_FINAL_ASSEMBLY);
1339: MatAssemblyEnd(ctx->Ybus,MAT_FINAL_ASSEMBLY);
1341: MatZeroEntries(ctx->J);
1343: ctx->alg_flg = PETSC_TRUE;
1345: /* Solve the algebraic equations */
1346: SNESSolve(snes_alg,NULL,X);
1348: ctx->stepnum++;
1350: /* Post-disturbance period */
1351: TSSetDuration(ts,1000,ctx->tmax);
1352: TSSetInitialTimeStep(ts,ctx->tfaultoff,.01);
1354: ctx->alg_flg = PETSC_TRUE;
1356: TSSolve(ts,X);
1357: TSGetTimeStepNumber(ts,&steps3);
1359: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1360: Adjoint model starts here
1361: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1362: TSSetPostStep(ts,NULL);
1363: MatCreateVecs(ctx->J,&lambda[0],NULL);
1364: /* Set initial conditions for the adjoint integration */
1365: VecZeroEntries(lambda[0]);
1367: MatCreateVecs(ctx->Jacp,&mu[0],NULL);
1368: VecZeroEntries(mu[0]);
1369: TSSetCostGradients(ts,1,lambda,mu);
1371: /* Set RHS JacobianP */
1372: TSAdjointSetRHSJacobian(ts,ctx->Jacp,RHSJacobianP,ctx);
1374: TSSetCostIntegrand(ts,1,(PetscErrorCode (*)(TS,PetscReal,Vec,Vec,void*))CostIntegrand,
1375: (PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*))DRDYFunction,
1376: (PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*))DRDPFunction,ctx);
1378: TSAdjointSetSteps(ts,steps3);
1379: TSAdjointSolve(ts);
1381: MatZeroEntries(ctx->J);
1382: /* Applying disturbance - resistive fault at ctx->faultbus */
1383: /* This is done by deducting shunt conductance to the diagonal location
1384: in the Ybus matrix */
1385: row_loc = 2*ctx->faultbus; col_loc = 2*ctx->faultbus+1; /* Location for G */
1386: val = 1./ctx->Rfault;
1387: MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);
1388: row_loc = 2*ctx->faultbus+1; col_loc = 2*ctx->faultbus; /* Location for G */
1389: val = 1./ctx->Rfault;
1390: MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);
1392: MatAssemblyBegin(ctx->Ybus,MAT_FINAL_ASSEMBLY);
1393: MatAssemblyEnd(ctx->Ybus,MAT_FINAL_ASSEMBLY);
1396: /* Set number of steps for the adjoint integration */
1397: TSAdjointSetSteps(ts,steps2);
1398: TSAdjointSolve(ts);
1400: MatZeroEntries(ctx->J);
1401: /* remove the fault */
1402: row_loc = 2*ctx->faultbus; col_loc = 2*ctx->faultbus+1; /* Location for G */
1403: val = -1./ctx->Rfault;
1404: MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);
1405: row_loc = 2*ctx->faultbus+1; col_loc = 2*ctx->faultbus; /* Location for G */
1406: val = -1./ctx->Rfault;
1407: MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);
1409: MatAssemblyBegin(ctx->Ybus,MAT_FINAL_ASSEMBLY);
1410: MatAssemblyEnd(ctx->Ybus,MAT_FINAL_ASSEMBLY);
1412: /* Set number of steps for the adjoint integration */
1413: TSAdjointSetSteps(ts,steps1);
1414: TSAdjointSolve(ts);
1417: ComputeSensiP(lambda[0],mu[0],DICDP,ctx);
1418: VecCopy(mu[0],G);
1419: TSGetCostIntegral(ts,&q);
1420: VecGetArray(q,&x_ptr);
1421: *f = x_ptr[0];
1423: VecDestroy(&lambda[0]);
1424: VecDestroy(&mu[0]);
1426: SNESDestroy(&snes_alg);
1427: VecDestroy(&F_alg);
1428: VecDestroy(&X);
1429: TSDestroy(&ts);
1430: for (i=0;i<3;i++) {
1431: VecDestroy(&DICDP[i]);
1432: }
1433: return 0;
1434: }