Actual source code: ex9busopt.c
petsc-3.7.1 2016-05-15
1: static char help[] = "Application of adjoint sensitivity analysis for 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: This code demonstrates how to solve a DAE-constrained optimization problem with TAO, TSAdjoint and TS.
10: The objectivie is to find optimal parameter PG for each generator to minizie the frequency violations due to faults.
11: The problem features discontinuities and a cost function in integral form.
12: The gradient is computed with the discrete adjoint of an implicit theta method, see ex9busadj.c for details.
13: */
15: #include <petsctao.h>
16: #include <petscts.h>
17: #include <petscdm.h>
18: #include <petscdmda.h>
19: #include <petscdmcomposite.h>
20: #include <petsctime.h>
22: PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal*,Vec,void*);
24: #define freq 60
25: #define w_s (2*PETSC_PI*freq)
27: /* Sizes and indices */
28: const PetscInt nbus = 9; /* Number of network buses */
29: const PetscInt ngen = 3; /* Number of generators */
30: const PetscInt nload = 3; /* Number of loads */
31: const PetscInt gbus[3] = {0,1,2}; /* Buses at which generators are incident */
32: const PetscInt lbus[3] = {4,5,7}; /* Buses at which loads are incident */
34: /* Generator real and reactive powers (found via loadflow) */
35: PetscScalar PG[3] = { 0.69,1.59,0.69};
36: /* PetscScalar PG[3] = {0.716786142395021,1.630000000000000,0.850000000000000};*/
38: const PetscScalar QG[3] = {0.270702180178785,0.066120127797275,-0.108402221791588};
39: /* Generator constants */
40: const PetscScalar H[3] = {23.64,6.4,3.01}; /* Inertia constant */
41: const PetscScalar Rs[3] = {0.0,0.0,0.0}; /* Stator Resistance */
42: const PetscScalar Xd[3] = {0.146,0.8958,1.3125}; /* d-axis reactance */
43: const PetscScalar Xdp[3] = {0.0608,0.1198,0.1813}; /* d-axis transient reactance */
44: 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 */
45: const PetscScalar Xqp[3] = {0.0969,0.1969,0.25}; /* q-axis transient reactance */
46: const PetscScalar Td0p[3] = {8.96,6.0,5.89}; /* d-axis open circuit time constant */
47: const PetscScalar Tq0p[3] = {0.31,0.535,0.6}; /* q-axis open circuit time constant */
48: PetscScalar M[3]; /* M = 2*H/w_s */
49: PetscScalar D[3]; /* D = 0.1*M */
51: PetscScalar TM[3]; /* Mechanical Torque */
52: /* Exciter system constants */
53: const PetscScalar KA[3] = {20.0,20.0,20.0}; /* Voltage regulartor gain constant */
54: const PetscScalar TA[3] = {0.2,0.2,0.2}; /* Voltage regulator time constant */
55: const PetscScalar KE[3] = {1.0,1.0,1.0}; /* Exciter gain constant */
56: const PetscScalar TE[3] = {0.314,0.314,0.314}; /* Exciter time constant */
57: const PetscScalar KF[3] = {0.063,0.063,0.063}; /* Feedback stabilizer gain constant */
58: const PetscScalar TF[3] = {0.35,0.35,0.35}; /* Feedback stabilizer time constant */
59: const PetscScalar k1[3] = {0.0039,0.0039,0.0039};
60: const PetscScalar k2[3] = {1.555,1.555,1.555}; /* k1 and k2 for calculating the saturation function SE = k1*exp(k2*Efd) */
62: PetscScalar Vref[3];
63: /* Load constants
64: We use a composite load model that describes the load and reactive powers at each time instant as follows
65: P(t) = \sum\limits_{i=0}^ld_nsegsp \ld_alphap_i*P_D0(\frac{V_m(t)}{V_m0})^\ld_betap_i
66: Q(t) = \sum\limits_{i=0}^ld_nsegsq \ld_alphaq_i*Q_D0(\frac{V_m(t)}{V_m0})^\ld_betaq_i
67: where
68: ld_nsegsp,ld_nsegsq - Number of individual load models for real and reactive power loads
69: ld_alphap,ld_alphap - Percentage contribution (weights) or loads
70: P_D0 - Real power load
71: Q_D0 - Reactive power load
72: V_m(t) - Voltage magnitude at time t
73: V_m0 - Voltage magnitude at t = 0
74: ld_betap, ld_betaq - exponents describing the load model for real and reactive part
76: Note: All loads have the same characteristic currently.
77: */
78: const PetscScalar PD0[3] = {1.25,0.9,1.0};
79: const PetscScalar QD0[3] = {0.5,0.3,0.35};
80: const PetscInt ld_nsegsp[3] = {3,3,3};
81: const PetscScalar ld_alphap[3] = {1.0,0.0,0.0};
82: const PetscScalar ld_betap[3] = {2.0,1.0,0.0};
83: const PetscInt ld_nsegsq[3] = {3,3,3};
84: const PetscScalar ld_alphaq[3] = {1.0,0.0,0.0};
85: const PetscScalar ld_betaq[3] = {2.0,1.0,0.0};
87: typedef struct {
88: DM dmgen, dmnet; /* DMs to manage generator and network subsystem */
89: DM dmpgrid; /* Composite DM to manage the entire power grid */
90: Mat Ybus; /* Network admittance matrix */
91: Vec V0; /* Initial voltage vector (Power flow solution) */
92: PetscReal tfaulton,tfaultoff; /* Fault on and off times */
93: PetscInt faultbus; /* Fault bus */
94: PetscScalar Rfault;
95: PetscReal t0,tmax;
96: PetscInt neqs_gen,neqs_net,neqs_pgrid;
97: Mat Sol; /* Matrix to save solution at each time step */
98: PetscInt stepnum;
99: PetscBool alg_flg;
100: PetscReal t;
101: IS is_diff; /* indices for differential equations */
102: IS is_alg; /* indices for algebraic equations */
103: PetscReal freq_u,freq_l; /* upper and lower frequency limit */
104: PetscInt pow; /* power coefficient used in the cost function */
105: PetscBool jacp_flg;
106: Mat J,Jacp;
107: } Userctx;
110: /* Converts from machine frame (dq) to network (phase a real,imag) reference frame */
113: PetscErrorCode dq2ri(PetscScalar Fd,PetscScalar Fq,PetscScalar delta,PetscScalar *Fr, PetscScalar *Fi)
114: {
116: *Fr = Fd*PetscSinScalar(delta) + Fq*PetscCosScalar(delta);
117: *Fi = -Fd*PetscCosScalar(delta) + Fq*PetscSinScalar(delta);
118: return(0);
119: }
121: /* Converts from network frame ([phase a real,imag) to machine (dq) reference frame */
124: PetscErrorCode ri2dq(PetscScalar Fr,PetscScalar Fi,PetscScalar delta,PetscScalar *Fd, PetscScalar *Fq)
125: {
127: *Fd = Fr*PetscSinScalar(delta) - Fi*PetscCosScalar(delta);
128: *Fq = Fr*PetscCosScalar(delta) + Fi*PetscSinScalar(delta);
129: return(0);
130: }
132: /* Saves the solution at each time to a matrix */
135: PetscErrorCode SaveSolution(TS ts)
136: {
137: PetscErrorCode ierr;
138: Userctx *user;
139: Vec X;
140: PetscScalar *mat;
141: const PetscScalar *x;
142: PetscInt idx;
143: PetscReal t;
146: TSGetApplicationContext(ts,&user);
147: TSGetTime(ts,&t);
148: TSGetSolution(ts,&X);
149: idx = user->stepnum*(user->neqs_pgrid+1);
150: MatDenseGetArray(user->Sol,&mat);
151: VecGetArrayRead(X,&x);
152: mat[idx] = t;
153: PetscMemcpy(mat+idx+1,x,user->neqs_pgrid*sizeof(PetscScalar));
154: MatDenseRestoreArray(user->Sol,&mat);
155: VecRestoreArrayRead(X,&x);
156: user->stepnum++;
157: return(0);
158: }
162: PetscErrorCode SetInitialGuess(Vec X,Userctx *user)
163: {
165: Vec Xgen,Xnet;
166: PetscScalar *xgen,*xnet;
167: PetscInt i,idx=0;
168: PetscScalar Vr,Vi,IGr,IGi,Vm,Vm2;
169: PetscScalar Eqp,Edp,delta;
170: PetscScalar Efd,RF,VR; /* Exciter variables */
171: PetscScalar Id,Iq; /* Generator dq axis currents */
172: PetscScalar theta,Vd,Vq,SE;
175: M[0] = 2*H[0]/w_s; M[1] = 2*H[1]/w_s; M[2] = 2*H[2]/w_s;
176: D[0] = 0.1*M[0]; D[1] = 0.1*M[1]; D[2] = 0.1*M[2];
178: DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
180: /* Network subsystem initialization */
181: VecCopy(user->V0,Xnet);
183: /* Generator subsystem initialization */
184: VecGetArray(Xgen,&xgen);
185: VecGetArray(Xnet,&xnet);
187: for (i=0; i < ngen; i++) {
188: Vr = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */
189: Vi = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */
190: Vm = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm;
191: IGr = (Vr*PG[i] + Vi*QG[i])/Vm2;
192: IGi = (Vi*PG[i] - Vr*QG[i])/Vm2;
194: delta = atan2(Vi+Xq[i]*IGr,Vr-Xq[i]*IGi); /* Machine angle */
196: theta = PETSC_PI/2.0 - delta;
198: Id = IGr*PetscCosScalar(theta) - IGi*PetscSinScalar(theta); /* d-axis stator current */
199: Iq = IGr*PetscSinScalar(theta) + IGi*PetscCosScalar(theta); /* q-axis stator current */
201: Vd = Vr*PetscCosScalar(theta) - Vi*PetscSinScalar(theta);
202: Vq = Vr*PetscSinScalar(theta) + Vi*PetscCosScalar(theta);
204: Edp = Vd + Rs[i]*Id - Xqp[i]*Iq; /* d-axis transient EMF */
205: Eqp = Vq + Rs[i]*Iq + Xdp[i]*Id; /* q-axis transient EMF */
207: TM[i] = PG[i];
209: /* The generator variables are ordered as [Eqp,Edp,delta,w,Id,Iq] */
210: xgen[idx] = Eqp;
211: xgen[idx+1] = Edp;
212: xgen[idx+2] = delta;
213: xgen[idx+3] = w_s;
215: idx = idx + 4;
217: xgen[idx] = Id;
218: xgen[idx+1] = Iq;
220: idx = idx + 2;
222: /* Exciter */
223: Efd = Eqp + (Xd[i] - Xdp[i])*Id;
224: SE = k1[i]*PetscExpScalar(k2[i]*Efd);
225: VR = KE[i]*Efd + SE;
226: RF = KF[i]*Efd/TF[i];
228: xgen[idx] = Efd;
229: xgen[idx+1] = RF;
230: xgen[idx+2] = VR;
232: Vref[i] = Vm + (VR/KA[i]);
234: idx = idx + 3;
235: }
237: VecRestoreArray(Xgen,&xgen);
238: VecRestoreArray(Xnet,&xnet);
240: /* VecView(Xgen,0); */
241: DMCompositeGather(user->dmpgrid,X,INSERT_VALUES,Xgen,Xnet);
242: DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);
243: return(0);
244: }
248: PetscErrorCode InitialGuess(Vec X,Userctx *user, const PetscScalar PGv[])
249: {
251: Vec Xgen,Xnet;
252: PetscScalar *xgen,*xnet;
253: PetscInt i,idx=0;
254: PetscScalar Vr,Vi,IGr,IGi,Vm,Vm2;
255: PetscScalar Eqp,Edp,delta;
256: PetscScalar Efd,RF,VR; /* Exciter variables */
257: PetscScalar Id,Iq; /* Generator dq axis currents */
258: PetscScalar theta,Vd,Vq,SE;
261: M[0] = 2*H[0]/w_s; M[1] = 2*H[1]/w_s; M[2] = 2*H[2]/w_s;
262: D[0] = 0.1*M[0]; D[1] = 0.1*M[1]; D[2] = 0.1*M[2];
264: DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
266: /* Network subsystem initialization */
267: VecCopy(user->V0,Xnet);
269: /* Generator subsystem initialization */
270: VecGetArray(Xgen,&xgen);
271: VecGetArray(Xnet,&xnet);
273: for (i=0; i < ngen; i++) {
274: Vr = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */
275: Vi = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */
276: Vm = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm;
277: IGr = (Vr*PGv[i] + Vi*QG[i])/Vm2;
278: IGi = (Vi*PGv[i] - Vr*QG[i])/Vm2;
280: delta = atan2(Vi+Xq[i]*IGr,Vr-Xq[i]*IGi); /* Machine angle */
282: theta = PETSC_PI/2.0 - delta;
284: Id = IGr*PetscCosScalar(theta) - IGi*PetscSinScalar(theta); /* d-axis stator current */
285: Iq = IGr*PetscSinScalar(theta) + IGi*PetscCosScalar(theta); /* q-axis stator current */
287: Vd = Vr*PetscCosScalar(theta) - Vi*PetscSinScalar(theta);
288: Vq = Vr*PetscSinScalar(theta) + Vi*PetscCosScalar(theta);
290: Edp = Vd + Rs[i]*Id - Xqp[i]*Iq; /* d-axis transient EMF */
291: Eqp = Vq + Rs[i]*Iq + Xdp[i]*Id; /* q-axis transient EMF */
293: /* The generator variables are ordered as [Eqp,Edp,delta,w,Id,Iq] */
294: xgen[idx] = Eqp;
295: xgen[idx+1] = Edp;
296: xgen[idx+2] = delta;
297: xgen[idx+3] = w_s;
299: idx = idx + 4;
301: xgen[idx] = Id;
302: xgen[idx+1] = Iq;
304: idx = idx + 2;
306: /* Exciter */
307: Efd = Eqp + (Xd[i] - Xdp[i])*Id;
308: SE = k1[i]*PetscExpScalar(k2[i]*Efd);
309: VR = KE[i]*Efd + SE;
310: RF = KF[i]*Efd/TF[i];
312: xgen[idx] = Efd;
313: xgen[idx+1] = RF;
314: xgen[idx+2] = VR;
316: idx = idx + 3;
317: }
319: VecRestoreArray(Xgen,&xgen);
320: VecRestoreArray(Xnet,&xnet);
322: /* VecView(Xgen,0); */
323: DMCompositeGather(user->dmpgrid,X,INSERT_VALUES,Xgen,Xnet);
324: DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);
325: return(0);
326: }
330: PetscErrorCode DICDPFiniteDifference(Vec X,Vec *DICDP, Userctx *user)
331: {
332: Vec Y;
333: PetscScalar PGv[3],eps;
335: PetscInt i,j;
337: eps = 1.e-7;
338: VecDuplicate(X,&Y);
340: for (i=0;i<ngen;i++) {
341: for (j=0;j<3;j++) PGv[j] = PG[j];
342: PGv[i] = PG[i]+eps;
343: InitialGuess(Y,user,PGv);
344: InitialGuess(X,user,PG);
346: VecAXPY(Y,-1.0,X);
347: VecScale(Y,1./eps);
348: VecCopy(Y,DICDP[i]);
349: }
350: VecDestroy(&Y);
351: return(0);
352: }
355: /* Computes F = [-f(x,y);g(x,y)] */
358: PetscErrorCode ResidualFunction(SNES snes,Vec X, Vec F, Userctx *user)
359: {
361: Vec Xgen,Xnet,Fgen,Fnet;
362: PetscScalar *xgen,*xnet,*fgen,*fnet;
363: PetscInt i,idx=0;
364: PetscScalar Vr,Vi,Vm,Vm2;
365: PetscScalar Eqp,Edp,delta,w; /* Generator variables */
366: PetscScalar Efd,RF,VR; /* Exciter variables */
367: PetscScalar Id,Iq; /* Generator dq axis currents */
368: PetscScalar Vd,Vq,SE;
369: PetscScalar IGr,IGi,IDr,IDi;
370: PetscScalar Zdq_inv[4],det;
371: PetscScalar PD,QD,Vm0,*v0;
372: PetscInt k;
375: VecZeroEntries(F);
376: DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
377: DMCompositeGetLocalVectors(user->dmpgrid,&Fgen,&Fnet);
378: DMCompositeScatter(user->dmpgrid,X,Xgen,Xnet);
379: DMCompositeScatter(user->dmpgrid,F,Fgen,Fnet);
381: /* Network current balance residual IG + Y*V + IL = 0. Only YV is added here.
382: The generator current injection, IG, and load current injection, ID are added later
383: */
384: /* Note that the values in Ybus are stored assuming the imaginary current balance
385: equation is ordered first followed by real current balance equation for each bus.
386: Thus imaginary current contribution goes in location 2*i, and
387: real current contribution in 2*i+1
388: */
389: MatMult(user->Ybus,Xnet,Fnet);
391: VecGetArray(Xgen,&xgen);
392: VecGetArray(Xnet,&xnet);
393: VecGetArray(Fgen,&fgen);
394: VecGetArray(Fnet,&fnet);
396: /* Generator subsystem */
397: for (i=0; i < ngen; i++) {
398: Eqp = xgen[idx];
399: Edp = xgen[idx+1];
400: delta = xgen[idx+2];
401: w = xgen[idx+3];
402: Id = xgen[idx+4];
403: Iq = xgen[idx+5];
404: Efd = xgen[idx+6];
405: RF = xgen[idx+7];
406: VR = xgen[idx+8];
408: /* Generator differential equations */
409: fgen[idx] = (Eqp + (Xd[i] - Xdp[i])*Id - Efd)/Td0p[i];
410: fgen[idx+1] = (Edp - (Xq[i] - Xqp[i])*Iq)/Tq0p[i];
411: fgen[idx+2] = -w + w_s;
412: fgen[idx+3] = (-TM[i] + Edp*Id + Eqp*Iq + (Xqp[i] - Xdp[i])*Id*Iq + D[i]*(w - w_s))/M[i];
414: Vr = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */
415: Vi = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */
417: ri2dq(Vr,Vi,delta,&Vd,&Vq);
418: /* Algebraic equations for stator currents */
419: det = Rs[i]*Rs[i] + Xdp[i]*Xqp[i];
421: Zdq_inv[0] = Rs[i]/det;
422: Zdq_inv[1] = Xqp[i]/det;
423: Zdq_inv[2] = -Xdp[i]/det;
424: Zdq_inv[3] = Rs[i]/det;
426: fgen[idx+4] = Zdq_inv[0]*(-Edp + Vd) + Zdq_inv[1]*(-Eqp + Vq) + Id;
427: fgen[idx+5] = Zdq_inv[2]*(-Edp + Vd) + Zdq_inv[3]*(-Eqp + Vq) + Iq;
429: /* Add generator current injection to network */
430: dq2ri(Id,Iq,delta,&IGr,&IGi);
432: fnet[2*gbus[i]] -= IGi;
433: fnet[2*gbus[i]+1] -= IGr;
435: Vm = PetscSqrtScalar(Vd*Vd + Vq*Vq); Vm2 = Vm*Vm;
437: SE = k1[i]*PetscExpScalar(k2[i]*Efd);
439: /* Exciter differential equations */
440: fgen[idx+6] = (KE[i]*Efd + SE - VR)/TE[i];
441: fgen[idx+7] = (RF - KF[i]*Efd/TF[i])/TF[i];
442: fgen[idx+8] = (VR - KA[i]*RF + KA[i]*KF[i]*Efd/TF[i] - KA[i]*(Vref[i] - Vm))/TA[i];
444: idx = idx + 9;
445: }
447: VecGetArray(user->V0,&v0);
448: for (i=0; i < nload; i++) {
449: Vr = xnet[2*lbus[i]]; /* Real part of load bus voltage */
450: Vi = xnet[2*lbus[i]+1]; /* Imaginary part of the load bus voltage */
451: Vm = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm;
452: Vm0 = PetscSqrtScalar(v0[2*lbus[i]]*v0[2*lbus[i]] + v0[2*lbus[i]+1]*v0[2*lbus[i]+1]);
453: PD = QD = 0.0;
454: for (k=0; k < ld_nsegsp[i]; k++) PD += ld_alphap[k]*PD0[i]*PetscPowScalar((Vm/Vm0),ld_betap[k]);
455: for (k=0; k < ld_nsegsq[i]; k++) QD += ld_alphaq[k]*QD0[i]*PetscPowScalar((Vm/Vm0),ld_betaq[k]);
457: /* Load currents */
458: IDr = (PD*Vr + QD*Vi)/Vm2;
459: IDi = (-QD*Vr + PD*Vi)/Vm2;
461: fnet[2*lbus[i]] += IDi;
462: fnet[2*lbus[i]+1] += IDr;
463: }
464: VecRestoreArray(user->V0,&v0);
466: VecRestoreArray(Xgen,&xgen);
467: VecRestoreArray(Xnet,&xnet);
468: VecRestoreArray(Fgen,&fgen);
469: VecRestoreArray(Fnet,&fnet);
471: DMCompositeGather(user->dmpgrid,F,INSERT_VALUES,Fgen,Fnet);
472: DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);
473: DMCompositeRestoreLocalVectors(user->dmpgrid,&Fgen,&Fnet);
474: return(0);
475: }
477: /* \dot{x} - f(x,y)
478: g(x,y) = 0
479: */
482: PetscErrorCode IFunction(TS ts,PetscReal t, Vec X, Vec Xdot, Vec F, Userctx *user)
483: {
484: PetscErrorCode ierr;
485: SNES snes;
486: PetscScalar *f;
487: const PetscScalar *xdot;
488: PetscInt i;
491: user->t = t;
493: TSGetSNES(ts,&snes);
494: ResidualFunction(snes,X,F,user);
495: VecGetArray(F,&f);
496: VecGetArrayRead(Xdot,&xdot);
497: for (i=0;i < ngen;i++) {
498: f[9*i] += xdot[9*i];
499: f[9*i+1] += xdot[9*i+1];
500: f[9*i+2] += xdot[9*i+2];
501: f[9*i+3] += xdot[9*i+3];
502: f[9*i+6] += xdot[9*i+6];
503: f[9*i+7] += xdot[9*i+7];
504: f[9*i+8] += xdot[9*i+8];
505: }
506: VecRestoreArray(F,&f);
507: VecRestoreArrayRead(Xdot,&xdot);
508: return(0);
509: }
511: /* This function is used for solving the algebraic system only during fault on and
512: off times. It computes the entire F and then zeros out the part corresponding to
513: differential equations
514: F = [0;g(y)];
515: */
518: PetscErrorCode AlgFunction(SNES snes, Vec X, Vec F, void *ctx)
519: {
521: Userctx *user=(Userctx*)ctx;
522: PetscScalar *f;
523: PetscInt i;
526: ResidualFunction(snes,X,F,user);
527: VecGetArray(F,&f);
528: for (i=0; i < ngen; i++) {
529: f[9*i] = 0;
530: f[9*i+1] = 0;
531: f[9*i+2] = 0;
532: f[9*i+3] = 0;
533: f[9*i+6] = 0;
534: f[9*i+7] = 0;
535: f[9*i+8] = 0;
536: }
537: VecRestoreArray(F,&f);
538: return(0);
539: }
543: PetscErrorCode PreallocateJacobian(Mat J, Userctx *user)
544: {
546: PetscInt *d_nnz;
547: PetscInt i,idx=0,start=0;
548: PetscInt ncols;
551: PetscMalloc1(user->neqs_pgrid,&d_nnz);
552: for (i=0; i<user->neqs_pgrid; i++) d_nnz[i] = 0;
553: /* Generator subsystem */
554: for (i=0; i < ngen; i++) {
556: d_nnz[idx] += 3;
557: d_nnz[idx+1] += 2;
558: d_nnz[idx+2] += 2;
559: d_nnz[idx+3] += 5;
560: d_nnz[idx+4] += 6;
561: d_nnz[idx+5] += 6;
563: d_nnz[user->neqs_gen+2*gbus[i]] += 3;
564: d_nnz[user->neqs_gen+2*gbus[i]+1] += 3;
566: d_nnz[idx+6] += 2;
567: d_nnz[idx+7] += 2;
568: d_nnz[idx+8] += 5;
570: idx = idx + 9;
571: }
573: start = user->neqs_gen;
574: for (i=0; i < nbus; i++) {
575: MatGetRow(user->Ybus,2*i,&ncols,NULL,NULL);
576: d_nnz[start+2*i] += ncols;
577: d_nnz[start+2*i+1] += ncols;
578: MatRestoreRow(user->Ybus,2*i,&ncols,NULL,NULL);
579: }
581: MatSeqAIJSetPreallocation(J,0,d_nnz);
582: PetscFree(d_nnz);
583: return(0);
584: }
586: /*
587: J = [-df_dx, -df_dy
588: dg_dx, dg_dy]
589: */
592: PetscErrorCode ResidualJacobian(SNES snes,Vec X,Mat J,Mat B,void *ctx)
593: {
594: PetscErrorCode ierr;
595: Userctx *user=(Userctx*)ctx;
596: Vec Xgen,Xnet;
597: PetscScalar *xgen,*xnet;
598: PetscInt i,idx=0;
599: PetscScalar Vr,Vi,Vm,Vm2;
600: PetscScalar Eqp,Edp,delta; /* Generator variables */
601: PetscScalar Efd; /* Exciter variables */
602: PetscScalar Id,Iq; /* Generator dq axis currents */
603: PetscScalar Vd,Vq;
604: PetscScalar val[10];
605: PetscInt row[2],col[10];
606: PetscInt net_start=user->neqs_gen;
607: PetscInt ncols;
608: const PetscInt *cols;
609: const PetscScalar *yvals;
610: PetscInt k;
611: PetscScalar Zdq_inv[4],det;
612: PetscScalar dVd_dVr,dVd_dVi,dVq_dVr,dVq_dVi,dVd_ddelta,dVq_ddelta;
613: PetscScalar dIGr_ddelta,dIGi_ddelta,dIGr_dId,dIGr_dIq,dIGi_dId,dIGi_dIq;
614: PetscScalar dSE_dEfd;
615: PetscScalar dVm_dVd,dVm_dVq,dVm_dVr,dVm_dVi;
616: PetscScalar PD,QD,Vm0,*v0,Vm4;
617: PetscScalar dPD_dVr,dPD_dVi,dQD_dVr,dQD_dVi;
618: PetscScalar dIDr_dVr,dIDr_dVi,dIDi_dVr,dIDi_dVi;
621: MatZeroEntries(B);
622: DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
623: DMCompositeScatter(user->dmpgrid,X,Xgen,Xnet);
625: VecGetArray(Xgen,&xgen);
626: VecGetArray(Xnet,&xnet);
628: /* Generator subsystem */
629: for (i=0; i < ngen; i++) {
630: Eqp = xgen[idx];
631: Edp = xgen[idx+1];
632: delta = xgen[idx+2];
633: Id = xgen[idx+4];
634: Iq = xgen[idx+5];
635: Efd = xgen[idx+6];
637: /* fgen[idx] = (Eqp + (Xd[i] - Xdp[i])*Id - Efd)/Td0p[i]; */
638: row[0] = idx;
639: col[0] = idx; col[1] = idx+4; col[2] = idx+6;
640: val[0] = 1/ Td0p[i]; val[1] = (Xd[i] - Xdp[i])/ Td0p[i]; val[2] = -1/Td0p[i];
642: MatSetValues(J,1,row,3,col,val,INSERT_VALUES);
644: /* fgen[idx+1] = (Edp - (Xq[i] - Xqp[i])*Iq)/Tq0p[i]; */
645: row[0] = idx + 1;
646: col[0] = idx + 1; col[1] = idx+5;
647: val[0] = 1/Tq0p[i]; val[1] = -(Xq[i] - Xqp[i])/Tq0p[i];
648: MatSetValues(J,1,row,2,col,val,INSERT_VALUES);
650: /* fgen[idx+2] = - w + w_s; */
651: row[0] = idx + 2;
652: col[0] = idx + 2; col[1] = idx + 3;
653: val[0] = 0; val[1] = -1;
654: MatSetValues(J,1,row,2,col,val,INSERT_VALUES);
656: /* fgen[idx+3] = (-TM[i] + Edp*Id + Eqp*Iq + (Xqp[i] - Xdp[i])*Id*Iq + D[i]*(w - w_s))/M[i]; */
657: row[0] = idx + 3;
658: col[0] = idx; col[1] = idx + 1; col[2] = idx + 3; col[3] = idx + 4; col[4] = idx + 5;
659: 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];
660: MatSetValues(J,1,row,5,col,val,INSERT_VALUES);
662: Vr = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */
663: Vi = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */
664: ri2dq(Vr,Vi,delta,&Vd,&Vq);
666: det = Rs[i]*Rs[i] + Xdp[i]*Xqp[i];
668: Zdq_inv[0] = Rs[i]/det;
669: Zdq_inv[1] = Xqp[i]/det;
670: Zdq_inv[2] = -Xdp[i]/det;
671: Zdq_inv[3] = Rs[i]/det;
673: dVd_dVr = PetscSinScalar(delta); dVd_dVi = -PetscCosScalar(delta);
674: dVq_dVr = PetscCosScalar(delta); dVq_dVi = PetscSinScalar(delta);
675: dVd_ddelta = Vr*PetscCosScalar(delta) + Vi*PetscSinScalar(delta);
676: dVq_ddelta = -Vr*PetscSinScalar(delta) + Vi*PetscCosScalar(delta);
678: /* fgen[idx+4] = Zdq_inv[0]*(-Edp + Vd) + Zdq_inv[1]*(-Eqp + Vq) + Id; */
679: row[0] = idx+4;
680: col[0] = idx; col[1] = idx+1; col[2] = idx + 2;
681: val[0] = -Zdq_inv[1]; val[1] = -Zdq_inv[0]; val[2] = Zdq_inv[0]*dVd_ddelta + Zdq_inv[1]*dVq_ddelta;
682: col[3] = idx + 4; col[4] = net_start+2*gbus[i]; col[5] = net_start + 2*gbus[i]+1;
683: 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;
684: MatSetValues(J,1,row,6,col,val,INSERT_VALUES);
686: /* fgen[idx+5] = Zdq_inv[2]*(-Edp + Vd) + Zdq_inv[3]*(-Eqp + Vq) + Iq; */
687: row[0] = idx+5;
688: col[0] = idx; col[1] = idx+1; col[2] = idx + 2;
689: val[0] = -Zdq_inv[3]; val[1] = -Zdq_inv[2]; val[2] = Zdq_inv[2]*dVd_ddelta + Zdq_inv[3]*dVq_ddelta;
690: col[3] = idx + 5; col[4] = net_start+2*gbus[i]; col[5] = net_start + 2*gbus[i]+1;
691: 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;
692: MatSetValues(J,1,row,6,col,val,INSERT_VALUES);
694: dIGr_ddelta = Id*PetscCosScalar(delta) - Iq*PetscSinScalar(delta);
695: dIGi_ddelta = Id*PetscSinScalar(delta) + Iq*PetscCosScalar(delta);
696: dIGr_dId = PetscSinScalar(delta); dIGr_dIq = PetscCosScalar(delta);
697: dIGi_dId = -PetscCosScalar(delta); dIGi_dIq = PetscSinScalar(delta);
699: /* fnet[2*gbus[i]] -= IGi; */
700: row[0] = net_start + 2*gbus[i];
701: col[0] = idx+2; col[1] = idx + 4; col[2] = idx + 5;
702: val[0] = -dIGi_ddelta; val[1] = -dIGi_dId; val[2] = -dIGi_dIq;
703: MatSetValues(J,1,row,3,col,val,INSERT_VALUES);
705: /* fnet[2*gbus[i]+1] -= IGr; */
706: row[0] = net_start + 2*gbus[i]+1;
707: col[0] = idx+2; col[1] = idx + 4; col[2] = idx + 5;
708: val[0] = -dIGr_ddelta; val[1] = -dIGr_dId; val[2] = -dIGr_dIq;
709: MatSetValues(J,1,row,3,col,val,INSERT_VALUES);
711: Vm = PetscSqrtScalar(Vd*Vd + Vq*Vq);
713: /* fgen[idx+6] = (KE[i]*Efd + SE - VR)/TE[i]; */
714: /* SE = k1[i]*PetscExpScalar(k2[i]*Efd); */
715: dSE_dEfd = k1[i]*k2[i]*PetscExpScalar(k2[i]*Efd);
717: row[0] = idx + 6;
718: col[0] = idx + 6; col[1] = idx + 8;
719: val[0] = (KE[i] + dSE_dEfd)/TE[i]; val[1] = -1/TE[i];
720: MatSetValues(J,1,row,2,col,val,INSERT_VALUES);
722: /* Exciter differential equations */
724: /* fgen[idx+7] = (RF - KF[i]*Efd/TF[i])/TF[i]; */
725: row[0] = idx + 7;
726: col[0] = idx + 6; col[1] = idx + 7;
727: val[0] = (-KF[i]/TF[i])/TF[i]; val[1] = 1/TF[i];
728: MatSetValues(J,1,row,2,col,val,INSERT_VALUES);
730: /* fgen[idx+8] = (VR - KA[i]*RF + KA[i]*KF[i]*Efd/TF[i] - KA[i]*(Vref[i] - Vm))/TA[i]; */
731: /* Vm = (Vd^2 + Vq^2)^0.5; */
732: dVm_dVd = Vd/Vm; dVm_dVq = Vq/Vm;
733: dVm_dVr = dVm_dVd*dVd_dVr + dVm_dVq*dVq_dVr;
734: dVm_dVi = dVm_dVd*dVd_dVi + dVm_dVq*dVq_dVi;
735: row[0] = idx + 8;
736: col[0] = idx + 6; col[1] = idx + 7; col[2] = idx + 8;
737: val[0] = (KA[i]*KF[i]/TF[i])/TA[i]; val[1] = -KA[i]/TA[i]; val[2] = 1/TA[i];
738: col[3] = net_start + 2*gbus[i]; col[4] = net_start + 2*gbus[i]+1;
739: val[3] = KA[i]*dVm_dVr/TA[i]; val[4] = KA[i]*dVm_dVi/TA[i];
740: MatSetValues(J,1,row,5,col,val,INSERT_VALUES);
741: idx = idx + 9;
742: }
745: for (i=0; i<nbus; i++) {
746: MatGetRow(user->Ybus,2*i,&ncols,&cols,&yvals);
747: row[0] = net_start + 2*i;
748: for (k=0; k<ncols; k++) {
749: col[k] = net_start + cols[k];
750: val[k] = yvals[k];
751: }
752: MatSetValues(J,1,row,ncols,col,val,INSERT_VALUES);
753: MatRestoreRow(user->Ybus,2*i,&ncols,&cols,&yvals);
755: MatGetRow(user->Ybus,2*i+1,&ncols,&cols,&yvals);
756: row[0] = net_start + 2*i+1;
757: for (k=0; k<ncols; k++) {
758: col[k] = net_start + cols[k];
759: val[k] = yvals[k];
760: }
761: MatSetValues(J,1,row,ncols,col,val,INSERT_VALUES);
762: MatRestoreRow(user->Ybus,2*i+1,&ncols,&cols,&yvals);
763: }
765: MatAssemblyBegin(J,MAT_FLUSH_ASSEMBLY);
766: MatAssemblyEnd(J,MAT_FLUSH_ASSEMBLY);
769: VecGetArray(user->V0,&v0);
770: for (i=0; i < nload; i++) {
771: Vr = xnet[2*lbus[i]]; /* Real part of load bus voltage */
772: Vi = xnet[2*lbus[i]+1]; /* Imaginary part of the load bus voltage */
773: Vm = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2= Vm*Vm; Vm4 = Vm2*Vm2;
774: Vm0 = PetscSqrtScalar(v0[2*lbus[i]]*v0[2*lbus[i]] + v0[2*lbus[i]+1]*v0[2*lbus[i]+1]);
775: PD = QD = 0.0;
776: dPD_dVr = dPD_dVi = dQD_dVr = dQD_dVi = 0.0;
777: for (k=0; k < ld_nsegsp[i]; k++) {
778: PD += ld_alphap[k]*PD0[i]*PetscPowScalar((Vm/Vm0),ld_betap[k]);
779: dPD_dVr += ld_alphap[k]*ld_betap[k]*PD0[i]*PetscPowScalar((1/Vm0),ld_betap[k])*Vr*PetscPowScalar(Vm,(ld_betap[k]-2));
780: dPD_dVi += ld_alphap[k]*ld_betap[k]*PD0[i]*PetscPowScalar((1/Vm0),ld_betap[k])*Vi*PetscPowScalar(Vm,(ld_betap[k]-2));
781: }
782: for (k=0; k < ld_nsegsq[i]; k++) {
783: QD += ld_alphaq[k]*QD0[i]*PetscPowScalar((Vm/Vm0),ld_betaq[k]);
784: dQD_dVr += ld_alphaq[k]*ld_betaq[k]*QD0[i]*PetscPowScalar((1/Vm0),ld_betaq[k])*Vr*PetscPowScalar(Vm,(ld_betaq[k]-2));
785: dQD_dVi += ld_alphaq[k]*ld_betaq[k]*QD0[i]*PetscPowScalar((1/Vm0),ld_betaq[k])*Vi*PetscPowScalar(Vm,(ld_betaq[k]-2));
786: }
788: /* IDr = (PD*Vr + QD*Vi)/Vm2; */
789: /* IDi = (-QD*Vr + PD*Vi)/Vm2; */
791: dIDr_dVr = (dPD_dVr*Vr + dQD_dVr*Vi + PD)/Vm2 - ((PD*Vr + QD*Vi)*2*Vr)/Vm4;
792: dIDr_dVi = (dPD_dVi*Vr + dQD_dVi*Vi + QD)/Vm2 - ((PD*Vr + QD*Vi)*2*Vi)/Vm4;
794: dIDi_dVr = (-dQD_dVr*Vr + dPD_dVr*Vi - QD)/Vm2 - ((-QD*Vr + PD*Vi)*2*Vr)/Vm4;
795: dIDi_dVi = (-dQD_dVi*Vr + dPD_dVi*Vi + PD)/Vm2 - ((-QD*Vr + PD*Vi)*2*Vi)/Vm4;
798: /* fnet[2*lbus[i]] += IDi; */
799: row[0] = net_start + 2*lbus[i];
800: col[0] = net_start + 2*lbus[i]; col[1] = net_start + 2*lbus[i]+1;
801: val[0] = dIDi_dVr; val[1] = dIDi_dVi;
802: MatSetValues(J,1,row,2,col,val,ADD_VALUES);
803: /* fnet[2*lbus[i]+1] += IDr; */
804: row[0] = net_start + 2*lbus[i]+1;
805: col[0] = net_start + 2*lbus[i]; col[1] = net_start + 2*lbus[i]+1;
806: val[0] = dIDr_dVr; val[1] = dIDr_dVi;
807: MatSetValues(J,1,row,2,col,val,ADD_VALUES);
808: }
809: VecRestoreArray(user->V0,&v0);
811: VecRestoreArray(Xgen,&xgen);
812: VecRestoreArray(Xnet,&xnet);
814: DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);
816: MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
817: MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
818: return(0);
819: }
821: /*
822: J = [I, 0
823: dg_dx, dg_dy]
824: */
827: PetscErrorCode AlgJacobian(SNES snes,Vec X,Mat A,Mat B,void *ctx)
828: {
830: Userctx *user=(Userctx*)ctx;
833: ResidualJacobian(snes,X,A,B,ctx);
834: MatSetOption(A,MAT_KEEP_NONZERO_PATTERN,PETSC_TRUE);
835: MatZeroRowsIS(A,user->is_diff,1.0,NULL,NULL);
836: return(0);
837: }
839: /*
840: J = [a*I-df_dx, -df_dy
841: dg_dx, dg_dy]
842: */
846: PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,Userctx *user)
847: {
849: SNES snes;
850: PetscScalar atmp = (PetscScalar) a;
851: PetscInt i,row;
854: user->t = t;
856: TSGetSNES(ts,&snes);
857: ResidualJacobian(snes,X,A,B,user);
858: for (i=0;i < ngen;i++) {
859: row = 9*i;
860: MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
861: row = 9*i+1;
862: MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
863: row = 9*i+2;
864: MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
865: row = 9*i+3;
866: MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
867: row = 9*i+6;
868: MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
869: row = 9*i+7;
870: MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
871: row = 9*i+8;
872: MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
873: }
874: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
875: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
876: return(0);
877: }
879: /* Matrix JacobianP is constant so that it only needs to be evaluated once */
882: static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A, void *ctx0)
883: {
885: PetscScalar a;
886: PetscInt row,col;
887: Userctx *ctx=(Userctx*)ctx0;
891: if (ctx->jacp_flg) {
892: MatZeroEntries(A);
894: for (col=0;col<3;col++) {
895: a = 1.0/M[col];
896: row = 9*col+3;
897: MatSetValues(A,1,&row,1,&col,&a,INSERT_VALUES);
898: }
900: ctx->jacp_flg = PETSC_FALSE;
902: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
903: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
904: }
905: return(0);
906: }
910: static PetscErrorCode CostIntegrand(TS ts,PetscReal t,Vec U,Vec R,Userctx *user)
911: {
913: PetscScalar *u,*r;
914: PetscInt idx;
915: Vec Xgen,Xnet;
916: PetscScalar *xgen;
917: PetscInt i;
920: DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
921: DMCompositeScatter(user->dmpgrid,U,Xgen,Xnet);
923: VecGetArray(Xgen,&xgen);
925: VecGetArray(U,&u);
926: VecGetArray(R,&r);
927: r[0] = 0.;
928: idx = 0;
929: for (i=0;i<ngen;i++) {
930: 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);
931: idx += 9;
932: }
933: VecRestoreArray(R,&r);
934: DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);
935: return(0);
936: }
940: static PetscErrorCode DRDYFunction(TS ts,PetscReal t,Vec U,Vec *drdy,Userctx *user)
941: {
943: Vec Xgen,Xnet,Dgen,Dnet;
944: PetscScalar *xgen,*dgen;
945: PetscInt i;
946: PetscInt idx;
949: DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
950: DMCompositeGetLocalVectors(user->dmpgrid,&Dgen,&Dnet);
951: DMCompositeScatter(user->dmpgrid,U,Xgen,Xnet);
952: DMCompositeScatter(user->dmpgrid,drdy[0],Dgen,Dnet);
954: VecGetArray(Xgen,&xgen);
955: VecGetArray(Dgen,&dgen);
957: idx = 0;
958: for (i=0;i<ngen;i++) {
959: dgen[idx+3] = 0.;
960: 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);
961: 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);
962: idx += 9;
963: }
965: VecRestoreArray(Dgen,&dgen);
966: DMCompositeGather(user->dmpgrid,drdy[0],INSERT_VALUES,Dgen,Dnet);
967: DMCompositeRestoreLocalVectors(user->dmpgrid,&Dgen,&Dnet);
968: DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);
969: return(0);
970: }
974: static PetscErrorCode DRDPFunction(TS ts,PetscReal t,Vec U,Vec *drdp,Userctx *user)
975: {
977: return(0);
978: }
982: PetscErrorCode ComputeSensiP(Vec lambda,Vec mu,Vec *DICDP,Userctx *user)
983: {
985: PetscScalar *x,*y,sensip;
986: PetscInt i;
989: VecGetArray(lambda,&x);
990: VecGetArray(mu,&y);
992: for (i=0;i<3;i++) {
993: VecDot(lambda,DICDP[i],&sensip);
994: sensip = sensip+y[i];
995: /* PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt %D th parameter: %g \n",i,(double)sensip); */
996: y[i] = sensip;
997: }
998: VecRestoreArray(mu,&y);
999: return(0);
1000: }
1004: int main(int argc,char **argv)
1005: {
1006: Userctx user;
1007: Vec p;
1008: PetscScalar *x_ptr;
1009: PetscErrorCode ierr;
1010: PetscMPIInt size;
1011: PetscInt i;
1012: PetscViewer Xview,Ybusview;
1013: PetscInt *idx2;
1014: Tao tao;
1015: KSP ksp;
1016: PC pc;
1017: Vec lowerb,upperb;
1019: PetscInitialize(&argc,&argv,"petscoptions",help);
1020: MPI_Comm_size(PETSC_COMM_WORLD,&size);
1021: if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs");
1023: user.jacp_flg = PETSC_TRUE;
1024: user.neqs_gen = 9*ngen; /* # eqs. for generator subsystem */
1025: user.neqs_net = 2*nbus; /* # eqs. for network subsystem */
1026: user.neqs_pgrid = user.neqs_gen + user.neqs_net;
1028: /* Create indices for differential and algebraic equations */
1029: PetscMalloc1(7*ngen,&idx2);
1030: for (i=0; i<ngen; i++) {
1031: 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;
1032: idx2[7*i+4] = 9*i+6; idx2[7*i+5] = 9*i+7; idx2[7*i+6] = 9*i+8;
1033: }
1034: ISCreateGeneral(PETSC_COMM_WORLD,7*ngen,idx2,PETSC_COPY_VALUES,&user.is_diff);
1035: ISComplement(user.is_diff,0,user.neqs_pgrid,&user.is_alg);
1036: PetscFree(idx2);
1038: /* Set run time options */
1039: PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Transient stability fault options","");
1040: {
1041: user.tfaulton = 1.0;
1042: user.tfaultoff = 1.2;
1043: user.Rfault = 0.0001;
1044: user.faultbus = 8;
1045: PetscOptionsReal("-tfaulton","","",user.tfaulton,&user.tfaulton,NULL);
1046: PetscOptionsReal("-tfaultoff","","",user.tfaultoff,&user.tfaultoff,NULL);
1047: PetscOptionsInt("-faultbus","","",user.faultbus,&user.faultbus,NULL);
1048: user.t0 = 0.0;
1049: user.tmax = 1.3;
1050: PetscOptionsReal("-t0","","",user.t0,&user.t0,NULL);
1051: PetscOptionsReal("-tmax","","",user.tmax,&user.tmax,NULL);
1052: user.freq_u = 61.0;
1053: user.freq_l = 59.0;
1054: user.pow = 2;
1055: PetscOptionsReal("-frequ","","",user.freq_u,&user.freq_u,NULL);
1056: PetscOptionsReal("-freql","","",user.freq_l,&user.freq_l,NULL);
1057: PetscOptionsInt("-pow","","",user.pow,&user.pow,NULL);
1059: }
1060: PetscOptionsEnd();
1062: /* Create DMs for generator and network subsystems */
1063: DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,user.neqs_gen,1,1,NULL,&user.dmgen);
1064: DMSetOptionsPrefix(user.dmgen,"dmgen_");
1065: DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,user.neqs_net,1,1,NULL,&user.dmnet);
1066: DMSetOptionsPrefix(user.dmnet,"dmnet_");
1067: /* Create a composite DM packer and add the two DMs */
1068: DMCompositeCreate(PETSC_COMM_WORLD,&user.dmpgrid);
1069: DMSetOptionsPrefix(user.dmpgrid,"pgrid_");
1070: DMCompositeAddDM(user.dmpgrid,user.dmgen);
1071: DMCompositeAddDM(user.dmpgrid,user.dmnet);
1073: /* Read initial voltage vector and Ybus */
1074: PetscViewerBinaryOpen(PETSC_COMM_WORLD,"X.bin",FILE_MODE_READ,&Xview);
1075: PetscViewerBinaryOpen(PETSC_COMM_WORLD,"Ybus.bin",FILE_MODE_READ,&Ybusview);
1077: VecCreate(PETSC_COMM_WORLD,&user.V0);
1078: VecSetSizes(user.V0,PETSC_DECIDE,user.neqs_net);
1079: VecLoad(user.V0,Xview);
1081: MatCreate(PETSC_COMM_WORLD,&user.Ybus);
1082: MatSetSizes(user.Ybus,PETSC_DECIDE,PETSC_DECIDE,user.neqs_net,user.neqs_net);
1083: MatSetType(user.Ybus,MATBAIJ);
1084: /* MatSetBlockSize(ctx->Ybus,2); */
1085: MatLoad(user.Ybus,Ybusview);
1087: PetscViewerDestroy(&Xview);
1088: PetscViewerDestroy(&Ybusview);
1090: /* Allocate space for Jacobians */
1091: MatCreate(PETSC_COMM_WORLD,&user.J);
1092: MatSetSizes(user.J,PETSC_DECIDE,PETSC_DECIDE,user.neqs_pgrid,user.neqs_pgrid);
1093: MatSetFromOptions(user.J);
1094: PreallocateJacobian(user.J,&user);
1096: MatCreate(PETSC_COMM_WORLD,&user.Jacp);
1097: MatSetSizes(user.Jacp,PETSC_DECIDE,PETSC_DECIDE,user.neqs_pgrid,3);
1098: MatSetFromOptions(user.Jacp);
1099: MatSetUp(user.Jacp);
1100: MatZeroEntries(user.Jacp); /* initialize to zeros */
1102: /* Create TAO solver and set desired solution method */
1103: TaoCreate(PETSC_COMM_WORLD,&tao);
1104: TaoSetType(tao,TAOBLMVM);
1105: /*
1106: Optimization starts
1107: */
1108: /* Set initial solution guess */
1109: VecCreateSeq(PETSC_COMM_WORLD,3,&p);
1110: VecGetArray(p,&x_ptr);
1111: x_ptr[0] = PG[0]; x_ptr[1] = PG[1]; x_ptr[2] = PG[2];
1112: VecRestoreArray(p,&x_ptr);
1114: TaoSetInitialVector(tao,p);
1115: /* Set routine for function and gradient evaluation */
1116: TaoSetObjectiveAndGradientRoutine(tao,FormFunctionGradient,&user);
1118: /* Set bounds for the optimization */
1119: VecDuplicate(p,&lowerb);
1120: VecDuplicate(p,&upperb);
1121: VecGetArray(lowerb,&x_ptr);
1122: x_ptr[0] = 0.5; x_ptr[1] = 0.5; x_ptr[2] = 0.5;
1123: VecRestoreArray(lowerb,&x_ptr);
1124: VecGetArray(upperb,&x_ptr);
1125: x_ptr[0] = 2.0; x_ptr[1] = 2.0; x_ptr[2] = 2.0;
1126: VecRestoreArray(upperb,&x_ptr);
1127: TaoSetVariableBounds(tao,lowerb,upperb);
1129: /* Check for any TAO command line options */
1130: TaoSetFromOptions(tao);
1131: TaoGetKSP(tao,&ksp);
1132: if (ksp) {
1133: KSPGetPC(ksp,&pc);
1134: PCSetType(pc,PCNONE);
1135: }
1137: /* SOLVE THE APPLICATION */
1138: TaoSolve(tao);
1140: VecView(p,PETSC_VIEWER_STDOUT_WORLD);
1141: /* Free TAO data structures */
1142: TaoDestroy(&tao);
1144: DMDestroy(&user.dmgen);
1145: DMDestroy(&user.dmnet);
1146: DMDestroy(&user.dmpgrid);
1147: ISDestroy(&user.is_diff);
1148: ISDestroy(&user.is_alg);
1150: MatDestroy(&user.J);
1151: MatDestroy(&user.Jacp);
1152: MatDestroy(&user.Ybus);
1153: /* MatDestroy(&user.Sol); */
1154: VecDestroy(&user.V0);
1155: VecDestroy(&p);
1156: VecDestroy(&lowerb);
1157: VecDestroy(&upperb);
1158: PetscFinalize();
1159: return(0);
1160: }
1162: /* ------------------------------------------------------------------ */
1165: /*
1166: FormFunction - Evaluates the function and corresponding gradient.
1168: Input Parameters:
1169: tao - the Tao context
1170: X - the input vector
1171: ptr - optional user-defined context, as set by TaoSetObjectiveAndGradientRoutine()
1173: Output Parameters:
1174: f - the newly evaluated function
1175: G - the newly evaluated gradient
1176: */
1177: PetscErrorCode FormFunctionGradient(Tao tao,Vec P,PetscReal *f,Vec G,void *ctx0)
1178: {
1179: TS ts;
1180: SNES snes_alg;
1182: Userctx *ctx = (Userctx*)ctx0;
1183: Vec X;
1184: PetscInt i;
1185: /* sensitivity context */
1186: PetscScalar *x_ptr;
1187: Vec lambda[1],q;
1188: Vec mu[1];
1189: PetscInt steps1,steps2,steps3;
1190: Vec DICDP[3];
1191: Vec F_alg;
1192: PetscInt row_loc,col_loc;
1193: PetscScalar val;
1194: Vec Xdot;
1196: VecGetArray(P,&x_ptr);
1197: PG[0] = x_ptr[0];
1198: PG[1] = x_ptr[1];
1199: PG[2] = x_ptr[2];
1200: VecRestoreArray(P,&x_ptr);
1202: ctx->stepnum = 0;
1204: DMCreateGlobalVector(ctx->dmpgrid,&X);
1206: /* Create matrix to save solutions at each time step */
1207: /* MatCreateSeqDense(PETSC_COMM_SELF,ctx->neqs_pgrid+1,1002,NULL,&ctx->Sol); */
1208: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1209: Create timestepping solver context
1210: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1211: TSCreate(PETSC_COMM_WORLD,&ts);
1212: TSSetProblemType(ts,TS_NONLINEAR);
1213: TSSetType(ts,TSCN);
1214: TSSetIFunction(ts,NULL,(TSIFunction) IFunction,ctx);
1215: TSSetIJacobian(ts,ctx->J,ctx->J,(TSIJacobian)IJacobian,ctx);
1216: TSSetApplicationContext(ts,ctx);
1218: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1219: Set initial conditions
1220: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1221: SetInitialGuess(X,ctx);
1223: /* Approximate DICDP with finite difference, we want to zero out network variables */
1224: for (i=0;i<3;i++) {
1225: VecDuplicate(X,&DICDP[i]);
1226: }
1227: DICDPFiniteDifference(X,DICDP,ctx);
1229: VecDuplicate(X,&F_alg);
1230: SNESCreate(PETSC_COMM_WORLD,&snes_alg);
1231: SNESSetFunction(snes_alg,F_alg,AlgFunction,ctx);
1232: MatZeroEntries(ctx->J);
1233: SNESSetJacobian(snes_alg,ctx->J,ctx->J,AlgJacobian,ctx);
1234: SNESSetOptionsPrefix(snes_alg,"alg_");
1235: SNESSetFromOptions(snes_alg);
1236: ctx->alg_flg = PETSC_TRUE;
1237: /* Solve the algebraic equations */
1238: SNESSolve(snes_alg,NULL,X);
1240: /* Just to set up the Jacobian structure */
1241: VecDuplicate(X,&Xdot);
1242: IJacobian(ts,0.0,X,Xdot,0.0,ctx->J,ctx->J,ctx);
1243: VecDestroy(&Xdot);
1245: ctx->stepnum++;
1247: /*
1248: Save trajectory of solution so that TSAdjointSolve() may be used
1249: */
1250: TSSetSaveTrajectory(ts);
1252: TSSetDuration(ts,1000,ctx->tfaulton);
1253: TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
1254: TSSetInitialTimeStep(ts,0.0,0.01);
1255: TSSetFromOptions(ts);
1256: /* TSSetPostStep(ts,SaveSolution); */
1258: ctx->alg_flg = PETSC_FALSE;
1259: /* Prefault period */
1260: TSSolve(ts,X);
1261: TSGetTimeStepNumber(ts,&steps1);
1263: /* Create the nonlinear solver for solving the algebraic system */
1264: /* Note that although the algebraic system needs to be solved only for
1265: Idq and V, we reuse the entire system including xgen. The xgen
1266: variables are held constant by setting their residuals to 0 and
1267: putting a 1 on the Jacobian diagonal for xgen rows
1268: */
1269: MatZeroEntries(ctx->J);
1271: /* Apply disturbance - resistive fault at ctx->faultbus */
1272: /* This is done by adding shunt conductance to the diagonal location
1273: in the Ybus matrix */
1274: row_loc = 2*ctx->faultbus; col_loc = 2*ctx->faultbus+1; /* Location for G */
1275: val = 1/ctx->Rfault;
1276: MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);
1277: row_loc = 2*ctx->faultbus+1; col_loc = 2*ctx->faultbus; /* Location for G */
1278: val = 1/ctx->Rfault;
1279: MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);
1281: MatAssemblyBegin(ctx->Ybus,MAT_FINAL_ASSEMBLY);
1282: MatAssemblyEnd(ctx->Ybus,MAT_FINAL_ASSEMBLY);
1284: ctx->alg_flg = PETSC_TRUE;
1285: /* Solve the algebraic equations */
1286: SNESSolve(snes_alg,NULL,X);
1288: ctx->stepnum++;
1290: /* Disturbance period */
1291: TSSetDuration(ts,1000,ctx->tfaultoff);
1292: TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
1293: TSSetInitialTimeStep(ts,ctx->tfaulton,.01);
1295: ctx->alg_flg = PETSC_FALSE;
1297: TSSolve(ts,X);
1298: TSGetTimeStepNumber(ts,&steps2);
1300: /* Remove the fault */
1301: row_loc = 2*ctx->faultbus; col_loc = 2*ctx->faultbus+1;
1302: val = -1/ctx->Rfault;
1303: MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);
1304: row_loc = 2*ctx->faultbus+1; col_loc = 2*ctx->faultbus;
1305: val = -1/ctx->Rfault;
1306: MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);
1308: MatAssemblyBegin(ctx->Ybus,MAT_FINAL_ASSEMBLY);
1309: MatAssemblyEnd(ctx->Ybus,MAT_FINAL_ASSEMBLY);
1311: MatZeroEntries(ctx->J);
1313: ctx->alg_flg = PETSC_TRUE;
1315: /* Solve the algebraic equations */
1316: SNESSolve(snes_alg,NULL,X);
1318: ctx->stepnum++;
1320: /* Post-disturbance period */
1321: TSSetDuration(ts,1000,ctx->tmax);
1322: TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
1323: TSSetInitialTimeStep(ts,ctx->tfaultoff,.01);
1325: ctx->alg_flg = PETSC_TRUE;
1327: TSSolve(ts,X);
1328: TSGetTimeStepNumber(ts,&steps3);
1330: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1331: Adjoint model starts here
1332: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1333: TSSetPostStep(ts,NULL);
1334: MatCreateVecs(ctx->J,&lambda[0],NULL);
1335: /* Set initial conditions for the adjoint integration */
1336: VecZeroEntries(lambda[0]);
1338: MatCreateVecs(ctx->Jacp,&mu[0],NULL);
1339: VecZeroEntries(mu[0]);
1340: TSSetCostGradients(ts,1,lambda,mu);
1342: /* Set RHS JacobianP */
1343: TSAdjointSetRHSJacobian(ts,ctx->Jacp,RHSJacobianP,ctx);
1345: TSSetCostIntegrand(ts,1,(PetscErrorCode (*)(TS,PetscReal,Vec,Vec,void*))CostIntegrand,
1346: (PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*))DRDYFunction,
1347: (PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*))DRDPFunction,PETSC_FALSE,ctx);
1349: TSAdjointSetSteps(ts,steps3);
1350: TSAdjointSolve(ts);
1352: MatZeroEntries(ctx->J);
1353: /* Applying disturbance - resistive fault at ctx->faultbus */
1354: /* This is done by deducting shunt conductance to the diagonal location
1355: in the Ybus matrix */
1356: row_loc = 2*ctx->faultbus; col_loc = 2*ctx->faultbus+1; /* Location for G */
1357: val = 1./ctx->Rfault;
1358: MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);
1359: row_loc = 2*ctx->faultbus+1; col_loc = 2*ctx->faultbus; /* Location for G */
1360: val = 1./ctx->Rfault;
1361: MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);
1363: MatAssemblyBegin(ctx->Ybus,MAT_FINAL_ASSEMBLY);
1364: MatAssemblyEnd(ctx->Ybus,MAT_FINAL_ASSEMBLY);
1367: /* Set number of steps for the adjoint integration */
1368: TSAdjointSetSteps(ts,steps2);
1369: TSAdjointSolve(ts);
1371: MatZeroEntries(ctx->J);
1372: /* remove the fault */
1373: row_loc = 2*ctx->faultbus; col_loc = 2*ctx->faultbus+1; /* Location for G */
1374: val = -1./ctx->Rfault;
1375: MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);
1376: row_loc = 2*ctx->faultbus+1; col_loc = 2*ctx->faultbus; /* Location for G */
1377: val = -1./ctx->Rfault;
1378: MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);
1380: MatAssemblyBegin(ctx->Ybus,MAT_FINAL_ASSEMBLY);
1381: MatAssemblyEnd(ctx->Ybus,MAT_FINAL_ASSEMBLY);
1383: /* Set number of steps for the adjoint integration */
1384: TSAdjointSetSteps(ts,steps1);
1385: TSAdjointSolve(ts);
1388: ComputeSensiP(lambda[0],mu[0],DICDP,ctx);
1389: VecCopy(mu[0],G);
1390: TSGetCostIntegral(ts,&q);
1391: VecGetArray(q,&x_ptr);
1392: *f = x_ptr[0];
1394: VecDestroy(&lambda[0]);
1395: VecDestroy(&mu[0]);
1397: SNESDestroy(&snes_alg);
1398: VecDestroy(&F_alg);
1399: VecDestroy(&X);
1400: TSDestroy(&ts);
1401: for (i=0;i<3;i++) {
1402: VecDestroy(&DICDP[i]);
1403: }
1404: return 0;
1405: }