Actual source code: ex9busadj.c
petsc-3.6.4 2016-04-12
2: static char help[] = "Power grid stability analysis of WECC 9 bus system.\n\
3: This example is based on the 9-bus (node) example given in the book Power\n\
4: Systems Dynamics and Stability (Chapter 7) by P. Sauer and M. A. Pai.\n\
5: The power grid in this example consists of 9 buses (nodes), 3 generators,\n\
6: 3 loads, and 9 transmission lines. The network equations are written\n\
7: in current balance form using rectangular coordiantes.\n\n";
9: /*
10: The equations for the stability analysis are described by the DAE
12: \dot{x} = f(x,y,t)
13: 0 = g(x,y,t)
15: where the generators are described by differential equations, while the algebraic
16: constraints define the network equations.
18: The generators are modeled with a 4th order differential equation describing the electrical
19: and mechanical dynamics. Each generator also has an exciter system modeled by 3rd order
20: diff. eqns. describing the exciter, voltage regulator, and the feedback stabilizer
21: mechanism.
23: The network equations are described by nodal current balance equations.
24: I(x,y) - Y*V = 0
26: where:
27: I(x,y) is the current injected from generators and loads.
28: Y is the admittance matrix, and
29: V is the voltage vector
30: */
32: /*
33: Include "petscts.h" so that we can use TS solvers. Note that this
34: file automatically includes:
35: petscsys.h - base PETSc routines petscvec.h - vectors
36: petscmat.h - matrices
37: petscis.h - index sets petscksp.h - Krylov subspace methods
38: petscviewer.h - viewers petscpc.h - preconditioners
39: petscksp.h - linear solvers
40: */
41: #include <petscts.h>
42: #include <petscdm.h>
43: #include <petscdmda.h>
44: #include <petscdmcomposite.h>
46: #define freq 60
47: #define w_s (2*PETSC_PI*freq)
49: /* Sizes and indices */
50: const PetscInt nbus = 9; /* Number of network buses */
51: const PetscInt ngen = 3; /* Number of generators */
52: const PetscInt nload = 3; /* Number of loads */
53: const PetscInt gbus[3] = {0,1,2}; /* Buses at which generators are incident */
54: const PetscInt lbus[3] = {4,5,7}; /* Buses at which loads are incident */
56: /* Generator real and reactive powers (found via loadflow) */
57: const PetscScalar PG[3] = {0.716786142395021,1.630000000000000,0.850000000000000};
58: const PetscScalar QG[3] = {0.270702180178785,0.066120127797275,-0.108402221791588};
59: /* Generator constants */
60: const PetscScalar H[3] = {23.64,6.4,3.01}; /* Inertia constant */
61: const PetscScalar Rs[3] = {0.0,0.0,0.0}; /* Stator Resistance */
62: const PetscScalar Xd[3] = {0.146,0.8958,1.3125}; /* d-axis reactance */
63: const PetscScalar Xdp[3] = {0.0608,0.1198,0.1813}; /* d-axis transient reactance */
64: const PetscScalar Xq[3] = {0.4360,0.8645,1.2578}; /* q-axis reactance Xq(1) set to 0.4360, value given in text 0.0969 */
65: const PetscScalar Xqp[3] = {0.0969,0.1969,0.25}; /* q-axis transient reactance */
66: const PetscScalar Td0p[3] = {8.96,6.0,5.89}; /* d-axis open circuit time constant */
67: const PetscScalar Tq0p[3] = {0.31,0.535,0.6}; /* q-axis open circuit time constant */
68: PetscScalar M[3]; /* M = 2*H/w_s */
69: PetscScalar D[3]; /* D = 0.1*M */
71: PetscScalar TM[3]; /* Mechanical Torque */
72: /* Exciter system constants */
73: const PetscScalar KA[3] = {20.0,20.0,20.0}; /* Voltage regulartor gain constant */
74: const PetscScalar TA[3] = {0.2,0.2,0.2}; /* Voltage regulator time constant */
75: const PetscScalar KE[3] = {1.0,1.0,1.0}; /* Exciter gain constant */
76: const PetscScalar TE[3] = {0.314,0.314,0.314}; /* Exciter time constant */
77: const PetscScalar KF[3] = {0.063,0.063,0.063}; /* Feedback stabilizer gain constant */
78: const PetscScalar TF[3] = {0.35,0.35,0.35}; /* Feedback stabilizer time constant */
79: const PetscScalar k1[3] = {0.0039,0.0039,0.0039};
80: const PetscScalar k2[3] = {1.555,1.555,1.555}; /* k1 and k2 for calculating the saturation function SE = k1*exp(k2*Efd) */
82: PetscScalar Vref[3];
83: /* Load constants
84: We use a composite load model that describes the load and reactive powers at each time instant as follows
85: P(t) = \sum\limits_{i=0}^ld_nsegsp \ld_alphap_i*P_D0(\frac{V_m(t)}{V_m0})^\ld_betap_i
86: Q(t) = \sum\limits_{i=0}^ld_nsegsq \ld_alphaq_i*Q_D0(\frac{V_m(t)}{V_m0})^\ld_betaq_i
87: where
88: ld_nsegsp,ld_nsegsq - Number of individual load models for real and reactive power loads
89: ld_alphap,ld_alphap - Percentage contribution (weights) or loads
90: P_D0 - Real power load
91: Q_D0 - Reactive power load
92: V_m(t) - Voltage magnitude at time t
93: V_m0 - Voltage magnitude at t = 0
94: ld_betap, ld_betaq - exponents describing the load model for real and reactive part
96: Note: All loads have the same characteristic currently.
97: */
98: const PetscScalar PD0[3] = {1.25,0.9,1.0};
99: const PetscScalar QD0[3] = {0.5,0.3,0.35};
100: const PetscInt ld_nsegsp[3] = {3,3,3};
101: const PetscScalar ld_alphap[3] = {1.0,0.0,0.0};
102: const PetscScalar ld_betap[3] = {2.0,1.0,0.0};
103: const PetscInt ld_nsegsq[3] = {3,3,3};
104: const PetscScalar ld_alphaq[3] = {1.0,0.0,0.0};
105: const PetscScalar ld_betaq[3] = {2.0,1.0,0.0};
107: typedef struct {
108: DM dmgen, dmnet; /* DMs to manage generator and network subsystem */
109: DM dmpgrid; /* Composite DM to manage the entire power grid */
110: Mat Ybus; /* Network admittance matrix */
111: Vec V0; /* Initial voltage vector (Power flow solution) */
112: PetscReal tfaulton,tfaultoff; /* Fault on and off times */
113: PetscInt faultbus; /* Fault bus */
114: PetscScalar Rfault;
115: PetscReal t0,tmax;
116: PetscInt neqs_gen,neqs_net,neqs_pgrid;
117: PetscBool alg_flg;
118: PetscReal t;
119: IS is_diff; /* indices for differential equations */
120: IS is_alg; /* indices for algebraic equations */
121: } Userctx;
124: /* Converts from machine frame (dq) to network (phase a real,imag) reference frame */
127: PetscErrorCode dq2ri(PetscScalar Fd,PetscScalar Fq,PetscScalar delta,PetscScalar *Fr, PetscScalar *Fi)
128: {
130: *Fr = Fd*PetscSinScalar(delta) + Fq*PetscCosScalar(delta);
131: *Fi = -Fd*PetscCosScalar(delta) + Fq*PetscSinScalar(delta);
132: return(0);
133: }
135: /* Converts from network frame ([phase a real,imag) to machine (dq) reference frame */
138: PetscErrorCode ri2dq(PetscScalar Fr,PetscScalar Fi,PetscScalar delta,PetscScalar *Fd, PetscScalar *Fq)
139: {
141: *Fd = Fr*PetscSinScalar(delta) - Fi*PetscCosScalar(delta);
142: *Fq = Fr*PetscCosScalar(delta) + Fi*PetscSinScalar(delta);
143: return(0);
144: }
148: PetscErrorCode SetInitialGuess(Vec X,Userctx *user)
149: {
151: Vec Xgen,Xnet;
152: PetscScalar *xgen,*xnet;
153: PetscInt i,idx=0;
154: PetscScalar Vr,Vi,IGr,IGi,Vm,Vm2;
155: PetscScalar Eqp,Edp,delta;
156: PetscScalar Efd,RF,VR; /* Exciter variables */
157: PetscScalar Id,Iq; /* Generator dq axis currents */
158: PetscScalar theta,Vd,Vq,SE;
161: M[0] = 2*H[0]/w_s; M[1] = 2*H[1]/w_s; M[2] = 2*H[2]/w_s;
162: D[0] = 0.1*M[0]; D[1] = 0.1*M[1]; D[2] = 0.1*M[2];
164: DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
166: /* Network subsystem initialization */
167: VecCopy(user->V0,Xnet);
169: /* Generator subsystem initialization */
170: VecGetArray(Xgen,&xgen);
171: VecGetArray(Xnet,&xnet);
173: for (i=0; i < ngen; i++) {
174: Vr = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */
175: Vi = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */
176: Vm = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm;
177: IGr = (Vr*PG[i] + Vi*QG[i])/Vm2;
178: IGi = (Vi*PG[i] - Vr*QG[i])/Vm2;
180: delta = atan2(Vi+Xq[i]*IGr,Vr-Xq[i]*IGi); /* Machine angle */
182: theta = PETSC_PI/2.0 - delta;
184: Id = IGr*PetscCosScalar(theta) - IGi*PetscSinScalar(theta); /* d-axis stator current */
185: Iq = IGr*PetscSinScalar(theta) + IGi*PetscCosScalar(theta); /* q-axis stator current */
187: Vd = Vr*PetscCosScalar(theta) - Vi*PetscSinScalar(theta);
188: Vq = Vr*PetscSinScalar(theta) + Vi*PetscCosScalar(theta);
190: Edp = Vd + Rs[i]*Id - Xqp[i]*Iq; /* d-axis transient EMF */
191: Eqp = Vq + Rs[i]*Iq + Xdp[i]*Id; /* q-axis transient EMF */
193: TM[i] = PG[i];
195: /* The generator variables are ordered as [Eqp,Edp,delta,w,Id,Iq] */
196: xgen[idx] = Eqp;
197: xgen[idx+1] = Edp;
198: xgen[idx+2] = delta;
199: xgen[idx+3] = w_s;
201: idx = idx + 4;
203: xgen[idx] = Id;
204: xgen[idx+1] = Iq;
206: idx = idx + 2;
208: /* Exciter */
209: Efd = Eqp + (Xd[i] - Xdp[i])*Id;
210: SE = k1[i]*PetscExpScalar(k2[i]*Efd);
211: VR = KE[i]*Efd + SE;
212: RF = KF[i]*Efd/TF[i];
214: xgen[idx] = Efd;
215: xgen[idx+1] = RF;
216: xgen[idx+2] = VR;
218: Vref[i] = Vm + (VR/KA[i]);
220: idx = idx + 3;
221: }
223: VecRestoreArray(Xgen,&xgen);
224: VecRestoreArray(Xnet,&xnet);
226: /* VecView(Xgen,0); */
227: DMCompositeGather(user->dmpgrid,X,INSERT_VALUES,Xgen,Xnet);
228: DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);
229: return(0);
230: }
232: /* Computes F = [f(x,y);g(x,y)] */
235: PetscErrorCode ResidualFunction(SNES snes,Vec X, Vec F, Userctx *user)
236: {
238: Vec Xgen,Xnet,Fgen,Fnet;
239: PetscScalar *xgen,*xnet,*fgen,*fnet;
240: PetscInt i,idx=0;
241: PetscScalar Vr,Vi,Vm,Vm2;
242: PetscScalar Eqp,Edp,delta,w; /* Generator variables */
243: PetscScalar Efd,RF,VR; /* Exciter variables */
244: PetscScalar Id,Iq; /* Generator dq axis currents */
245: PetscScalar Vd,Vq,SE;
246: PetscScalar IGr,IGi,IDr,IDi;
247: PetscScalar Zdq_inv[4],det;
248: PetscScalar PD,QD,Vm0,*v0;
249: PetscInt k;
252: VecZeroEntries(F);
253: DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
254: DMCompositeGetLocalVectors(user->dmpgrid,&Fgen,&Fnet);
255: DMCompositeScatter(user->dmpgrid,X,Xgen,Xnet);
256: DMCompositeScatter(user->dmpgrid,F,Fgen,Fnet);
258: /* Network current balance residual IG + Y*V + IL = 0. Only YV is added here.
259: The generator current injection, IG, and load current injection, ID are added later
260: */
261: /* Note that the values in Ybus are stored assuming the imaginary current balance
262: equation is ordered first followed by real current balance equation for each bus.
263: Thus imaginary current contribution goes in location 2*i, and
264: real current contribution in 2*i+1
265: */
266: MatMult(user->Ybus,Xnet,Fnet);
268: VecGetArray(Xgen,&xgen);
269: VecGetArray(Xnet,&xnet);
270: VecGetArray(Fgen,&fgen);
271: VecGetArray(Fnet,&fnet);
273: /* Generator subsystem */
274: for (i=0; i < ngen; i++) {
275: Eqp = xgen[idx];
276: Edp = xgen[idx+1];
277: delta = xgen[idx+2];
278: w = xgen[idx+3];
279: Id = xgen[idx+4];
280: Iq = xgen[idx+5];
281: Efd = xgen[idx+6];
282: RF = xgen[idx+7];
283: VR = xgen[idx+8];
285: /* Generator differential equations */
286: fgen[idx] = (Eqp + (Xd[i] - Xdp[i])*Id - Efd)/Td0p[i];
287: fgen[idx+1] = (Edp - (Xq[i] - Xqp[i])*Iq)/Tq0p[i];
288: fgen[idx+2] = -w + w_s;
289: fgen[idx+3] = (-TM[i] + Edp*Id + Eqp*Iq + (Xqp[i] - Xdp[i])*Id*Iq + D[i]*(w - w_s))/M[i];
291: Vr = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */
292: Vi = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */
294: ri2dq(Vr,Vi,delta,&Vd,&Vq);
295: /* Algebraic equations for stator currents */
297: det = Rs[i]*Rs[i] + Xdp[i]*Xqp[i];
299: Zdq_inv[0] = Rs[i]/det;
300: Zdq_inv[1] = Xqp[i]/det;
301: Zdq_inv[2] = -Xdp[i]/det;
302: Zdq_inv[3] = Rs[i]/det;
304: fgen[idx+4] = Zdq_inv[0]*(-Edp + Vd) + Zdq_inv[1]*(-Eqp + Vq) + Id;
305: fgen[idx+5] = Zdq_inv[2]*(-Edp + Vd) + Zdq_inv[3]*(-Eqp + Vq) + Iq;
307: /* Add generator current injection to network */
308: dq2ri(Id,Iq,delta,&IGr,&IGi);
310: fnet[2*gbus[i]] -= IGi;
311: fnet[2*gbus[i]+1] -= IGr;
313: Vm = PetscSqrtScalar(Vd*Vd + Vq*Vq); Vm2 = Vm*Vm;
315: SE = k1[i]*PetscExpScalar(k2[i]*Efd);
317: /* Exciter differential equations */
318: fgen[idx+6] = (KE[i]*Efd + SE - VR)/TE[i];
319: fgen[idx+7] = (RF - KF[i]*Efd/TF[i])/TF[i];
320: fgen[idx+8] = (VR - KA[i]*RF + KA[i]*KF[i]*Efd/TF[i] - KA[i]*(Vref[i] - Vm))/TA[i];
322: idx = idx + 9;
323: }
325: VecGetArray(user->V0,&v0);
326: for (i=0; i < nload; i++) {
327: Vr = xnet[2*lbus[i]]; /* Real part of load bus voltage */
328: Vi = xnet[2*lbus[i]+1]; /* Imaginary part of the load bus voltage */
329: Vm = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm;
330: Vm0 = PetscSqrtScalar(v0[2*lbus[i]]*v0[2*lbus[i]] + v0[2*lbus[i]+1]*v0[2*lbus[i]+1]);
331: PD = QD = 0.0;
332: for (k=0; k < ld_nsegsp[i]; k++) PD += ld_alphap[k]*PD0[i]*PetscPowScalar((Vm/Vm0),ld_betap[k]);
333: for (k=0; k < ld_nsegsq[i]; k++) QD += ld_alphaq[k]*QD0[i]*PetscPowScalar((Vm/Vm0),ld_betaq[k]);
335: /* Load currents */
336: IDr = (PD*Vr + QD*Vi)/Vm2;
337: IDi = (-QD*Vr + PD*Vi)/Vm2;
339: fnet[2*lbus[i]] += IDi;
340: fnet[2*lbus[i]+1] += IDr;
341: }
342: VecRestoreArray(user->V0,&v0);
344: VecRestoreArray(Xgen,&xgen);
345: VecRestoreArray(Xnet,&xnet);
346: VecRestoreArray(Fgen,&fgen);
347: VecRestoreArray(Fnet,&fnet);
349: DMCompositeGather(user->dmpgrid,F,INSERT_VALUES,Fgen,Fnet);
350: DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);
351: DMCompositeRestoreLocalVectors(user->dmpgrid,&Fgen,&Fnet);
352: return(0);
353: }
355: /* \dot{x} - f(x,y)
356: g(x,y) = 0
357: */
360: PetscErrorCode IFunction(TS ts,PetscReal t, Vec X, Vec Xdot, Vec F, Userctx *user)
361: {
362: PetscErrorCode ierr;
363: SNES snes;
364: PetscScalar *f;
365: const PetscScalar *xdot;
366: PetscInt i;
369: user->t = t;
371: TSGetSNES(ts,&snes);
372: ResidualFunction(snes,X,F,user);
373: VecGetArray(F,&f);
374: VecGetArrayRead(Xdot,&xdot);
375: for (i=0;i < ngen;i++) {
376: f[9*i] += xdot[9*i];
377: f[9*i+1] += xdot[9*i+1];
378: f[9*i+2] += xdot[9*i+2];
379: f[9*i+3] += xdot[9*i+3];
380: f[9*i+6] += xdot[9*i+6];
381: f[9*i+7] += xdot[9*i+7];
382: f[9*i+8] += xdot[9*i+8];
383: }
384: VecRestoreArray(F,&f);
385: VecRestoreArrayRead(Xdot,&xdot);
386: return(0);
387: }
389: /* This function is used for solving the algebraic system only during fault on and
390: off times. It computes the entire F and then zeros out the part corresponding to
391: differential equations
392: F = [0;g(y)];
393: */
396: PetscErrorCode AlgFunction(SNES snes, Vec X, Vec F, void *ctx)
397: {
399: Userctx *user=(Userctx*)ctx;
400: PetscScalar *f;
401: PetscInt i;
404: ResidualFunction(snes,X,F,user);
405: VecGetArray(F,&f);
406: for (i=0; i < ngen; i++) {
407: f[9*i] = 0;
408: f[9*i+1] = 0;
409: f[9*i+2] = 0;
410: f[9*i+3] = 0;
411: f[9*i+6] = 0;
412: f[9*i+7] = 0;
413: f[9*i+8] = 0;
414: }
415: VecRestoreArray(F,&f);
416: return(0);
417: }
421: PetscErrorCode PreallocateJacobian(Mat J, Userctx *user)
422: {
424: PetscInt *d_nnz;
425: PetscInt i,idx=0,start=0;
426: PetscInt ncols;
429: PetscMalloc1(user->neqs_pgrid,&d_nnz);
430: for (i=0; i<user->neqs_pgrid; i++) d_nnz[i] = 0;
431: /* Generator subsystem */
432: for (i=0; i < ngen; i++) {
434: d_nnz[idx] += 3;
435: d_nnz[idx+1] += 2;
436: d_nnz[idx+2] += 2;
437: d_nnz[idx+3] += 5;
438: d_nnz[idx+4] += 6;
439: d_nnz[idx+5] += 6;
441: d_nnz[user->neqs_gen+2*gbus[i]] += 3;
442: d_nnz[user->neqs_gen+2*gbus[i]+1] += 3;
444: d_nnz[idx+6] += 2;
445: d_nnz[idx+7] += 2;
446: d_nnz[idx+8] += 5;
448: idx = idx + 9;
449: }
451: start = user->neqs_gen;
453: for (i=0; i < nbus; i++) {
454: MatGetRow(user->Ybus,2*i,&ncols,NULL,NULL);
455: d_nnz[start+2*i] += ncols;
456: d_nnz[start+2*i+1] += ncols;
457: MatRestoreRow(user->Ybus,2*i,&ncols,NULL,NULL);
458: }
460: MatSeqAIJSetPreallocation(J,0,d_nnz);
462: PetscFree(d_nnz);
463: return(0);
464: }
466: /*
467: J = [-df_dx, -df_dy
468: dg_dx, dg_dy]
469: */
472: PetscErrorCode ResidualJacobian(SNES snes,Vec X,Mat J,Mat B,void *ctx)
473: {
474: PetscErrorCode ierr;
475: Userctx *user=(Userctx*)ctx;
476: Vec Xgen,Xnet;
477: PetscScalar *xgen,*xnet;
478: PetscInt i,idx=0;
479: PetscScalar Vr,Vi,Vm,Vm2;
480: PetscScalar Eqp,Edp,delta; /* Generator variables */
481: PetscScalar Efd; /* Exciter variables */
482: PetscScalar Id,Iq; /* Generator dq axis currents */
483: PetscScalar Vd,Vq;
484: PetscScalar val[10];
485: PetscInt row[2],col[10];
486: PetscInt net_start=user->neqs_gen;
487: PetscScalar Zdq_inv[4],det;
488: PetscScalar dVd_dVr,dVd_dVi,dVq_dVr,dVq_dVi,dVd_ddelta,dVq_ddelta;
489: PetscScalar dIGr_ddelta,dIGi_ddelta,dIGr_dId,dIGr_dIq,dIGi_dId,dIGi_dIq;
490: PetscScalar dSE_dEfd;
491: PetscScalar dVm_dVd,dVm_dVq,dVm_dVr,dVm_dVi;
492: PetscInt ncols;
493: const PetscInt *cols;
494: const PetscScalar *yvals;
495: PetscInt k;
496: PetscScalar PD,QD,Vm0,*v0,Vm4;
497: PetscScalar dPD_dVr,dPD_dVi,dQD_dVr,dQD_dVi;
498: PetscScalar dIDr_dVr,dIDr_dVi,dIDi_dVr,dIDi_dVi;
501: MatZeroEntries(B);
502: DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
503: DMCompositeScatter(user->dmpgrid,X,Xgen,Xnet);
505: VecGetArray(Xgen,&xgen);
506: VecGetArray(Xnet,&xnet);
508: /* Generator subsystem */
509: for (i=0; i < ngen; i++) {
510: Eqp = xgen[idx];
511: Edp = xgen[idx+1];
512: delta = xgen[idx+2];
513: Id = xgen[idx+4];
514: Iq = xgen[idx+5];
515: Efd = xgen[idx+6];
517: /* fgen[idx] = (Eqp + (Xd[i] - Xdp[i])*Id - Efd)/Td0p[i]; */
518: row[0] = idx;
519: col[0] = idx; col[1] = idx+4; col[2] = idx+6;
520: val[0] = 1/ Td0p[i]; val[1] = (Xd[i] - Xdp[i])/ Td0p[i]; val[2] = -1/Td0p[i];
522: MatSetValues(J,1,row,3,col,val,INSERT_VALUES);
524: /* fgen[idx+1] = (Edp - (Xq[i] - Xqp[i])*Iq)/Tq0p[i]; */
525: row[0] = idx + 1;
526: col[0] = idx + 1; col[1] = idx+5;
527: val[0] = 1/Tq0p[i]; val[1] = -(Xq[i] - Xqp[i])/Tq0p[i];
528: MatSetValues(J,1,row,2,col,val,INSERT_VALUES);
530: /* fgen[idx+2] = - w + w_s; */
531: row[0] = idx + 2;
532: col[0] = idx + 2; col[1] = idx + 3;
533: val[0] = 0; val[1] = -1;
534: MatSetValues(J,1,row,2,col,val,INSERT_VALUES);
536: /* fgen[idx+3] = (-TM[i] + Edp*Id + Eqp*Iq + (Xqp[i] - Xdp[i])*Id*Iq + D[i]*(w - w_s))/M[i]; */
537: row[0] = idx + 3;
538: col[0] = idx; col[1] = idx + 1; col[2] = idx + 3; col[3] = idx + 4; col[4] = idx + 5;
539: 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];
540: MatSetValues(J,1,row,5,col,val,INSERT_VALUES);
542: Vr = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */
543: Vi = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */
544: ri2dq(Vr,Vi,delta,&Vd,&Vq);
546: det = Rs[i]*Rs[i] + Xdp[i]*Xqp[i];
548: Zdq_inv[0] = Rs[i]/det;
549: Zdq_inv[1] = Xqp[i]/det;
550: Zdq_inv[2] = -Xdp[i]/det;
551: Zdq_inv[3] = Rs[i]/det;
553: dVd_dVr = PetscSinScalar(delta); dVd_dVi = -PetscCosScalar(delta);
554: dVq_dVr = PetscCosScalar(delta); dVq_dVi = PetscSinScalar(delta);
555: dVd_ddelta = Vr*PetscCosScalar(delta) + Vi*PetscSinScalar(delta);
556: dVq_ddelta = -Vr*PetscSinScalar(delta) + Vi*PetscCosScalar(delta);
558: /* fgen[idx+4] = Zdq_inv[0]*(-Edp + Vd) + Zdq_inv[1]*(-Eqp + Vq) + Id; */
559: row[0] = idx+4;
560: col[0] = idx; col[1] = idx+1; col[2] = idx + 2;
561: val[0] = -Zdq_inv[1]; val[1] = -Zdq_inv[0]; val[2] = Zdq_inv[0]*dVd_ddelta + Zdq_inv[1]*dVq_ddelta;
562: col[3] = idx + 4; col[4] = net_start+2*gbus[i]; col[5] = net_start + 2*gbus[i]+1;
563: 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;
564: MatSetValues(J,1,row,6,col,val,INSERT_VALUES);
566: /* fgen[idx+5] = Zdq_inv[2]*(-Edp + Vd) + Zdq_inv[3]*(-Eqp + Vq) + Iq; */
567: row[0] = idx+5;
568: col[0] = idx; col[1] = idx+1; col[2] = idx + 2;
569: val[0] = -Zdq_inv[3]; val[1] = -Zdq_inv[2]; val[2] = Zdq_inv[2]*dVd_ddelta + Zdq_inv[3]*dVq_ddelta;
570: col[3] = idx + 5; col[4] = net_start+2*gbus[i]; col[5] = net_start + 2*gbus[i]+1;
571: 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;
572: MatSetValues(J,1,row,6,col,val,INSERT_VALUES);
574: dIGr_ddelta = Id*PetscCosScalar(delta) - Iq*PetscSinScalar(delta);
575: dIGi_ddelta = Id*PetscSinScalar(delta) + Iq*PetscCosScalar(delta);
576: dIGr_dId = PetscSinScalar(delta); dIGr_dIq = PetscCosScalar(delta);
577: dIGi_dId = -PetscCosScalar(delta); dIGi_dIq = PetscSinScalar(delta);
579: /* fnet[2*gbus[i]] -= IGi; */
580: row[0] = net_start + 2*gbus[i];
581: col[0] = idx+2; col[1] = idx + 4; col[2] = idx + 5;
582: val[0] = -dIGi_ddelta; val[1] = -dIGi_dId; val[2] = -dIGi_dIq;
583: MatSetValues(J,1,row,3,col,val,INSERT_VALUES);
585: /* fnet[2*gbus[i]+1] -= IGr; */
586: row[0] = net_start + 2*gbus[i]+1;
587: col[0] = idx+2; col[1] = idx + 4; col[2] = idx + 5;
588: val[0] = -dIGr_ddelta; val[1] = -dIGr_dId; val[2] = -dIGr_dIq;
589: MatSetValues(J,1,row,3,col,val,INSERT_VALUES);
591: Vm = PetscSqrtScalar(Vd*Vd + Vq*Vq); Vm2 = Vm*Vm;
593: /* fgen[idx+6] = (KE[i]*Efd + SE - VR)/TE[i]; */
594: /* SE = k1[i]*PetscExpScalar(k2[i]*Efd); */
595: dSE_dEfd = k1[i]*k2[i]*PetscExpScalar(k2[i]*Efd);
597: row[0] = idx + 6;
598: col[0] = idx + 6; col[1] = idx + 8;
599: val[0] = (KE[i] + dSE_dEfd)/TE[i]; val[1] = -1/TE[i];
600: MatSetValues(J,1,row,2,col,val,INSERT_VALUES);
602: /* Exciter differential equations */
604: /* fgen[idx+7] = (RF - KF[i]*Efd/TF[i])/TF[i]; */
605: row[0] = idx + 7;
606: col[0] = idx + 6; col[1] = idx + 7;
607: val[0] = (-KF[i]/TF[i])/TF[i]; val[1] = 1/TF[i];
608: MatSetValues(J,1,row,2,col,val,INSERT_VALUES);
610: /* fgen[idx+8] = (VR - KA[i]*RF + KA[i]*KF[i]*Efd/TF[i] - KA[i]*(Vref[i] - Vm))/TA[i]; */
611: /* Vm = (Vd^2 + Vq^2)^0.5; */
612: dVm_dVd = Vd/Vm; dVm_dVq = Vq/Vm;
613: dVm_dVr = dVm_dVd*dVd_dVr + dVm_dVq*dVq_dVr;
614: dVm_dVi = dVm_dVd*dVd_dVi + dVm_dVq*dVq_dVi;
615: row[0] = idx + 8;
616: col[0] = idx + 6; col[1] = idx + 7; col[2] = idx + 8;
617: val[0] = (KA[i]*KF[i]/TF[i])/TA[i]; val[1] = -KA[i]/TA[i]; val[2] = 1/TA[i];
618: col[3] = net_start + 2*gbus[i]; col[4] = net_start + 2*gbus[i]+1;
619: val[3] = KA[i]*dVm_dVr/TA[i]; val[4] = KA[i]*dVm_dVi/TA[i];
620: MatSetValues(J,1,row,5,col,val,INSERT_VALUES);
621: idx = idx + 9;
622: }
624: for (i=0; i<nbus; i++) {
625: MatGetRow(user->Ybus,2*i,&ncols,&cols,&yvals);
626: row[0] = net_start + 2*i;
627: for (k=0; k<ncols; k++) {
628: col[k] = net_start + cols[k];
629: val[k] = yvals[k];
630: }
631: MatSetValues(J,1,row,ncols,col,val,INSERT_VALUES);
632: MatRestoreRow(user->Ybus,2*i,&ncols,&cols,&yvals);
634: MatGetRow(user->Ybus,2*i+1,&ncols,&cols,&yvals);
635: row[0] = net_start + 2*i+1;
636: for (k=0; k<ncols; k++) {
637: col[k] = net_start + cols[k];
638: val[k] = yvals[k];
639: }
640: MatSetValues(J,1,row,ncols,col,val,INSERT_VALUES);
641: MatRestoreRow(user->Ybus,2*i+1,&ncols,&cols,&yvals);
642: }
644: MatAssemblyBegin(J,MAT_FLUSH_ASSEMBLY);
645: MatAssemblyEnd(J,MAT_FLUSH_ASSEMBLY);
648: VecGetArray(user->V0,&v0);
649: for (i=0; i < nload; i++) {
650: Vr = xnet[2*lbus[i]]; /* Real part of load bus voltage */
651: Vi = xnet[2*lbus[i]+1]; /* Imaginary part of the load bus voltage */
652: Vm = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm; Vm4 = Vm2*Vm2;
653: Vm0 = PetscSqrtScalar(v0[2*lbus[i]]*v0[2*lbus[i]] + v0[2*lbus[i]+1]*v0[2*lbus[i]+1]);
654: PD = QD = 0.0;
655: dPD_dVr = dPD_dVi = dQD_dVr = dQD_dVi = 0.0;
656: for (k=0; k < ld_nsegsp[i]; k++) {
657: PD += ld_alphap[k]*PD0[i]*PetscPowScalar((Vm/Vm0),ld_betap[k]);
658: dPD_dVr += ld_alphap[k]*ld_betap[k]*PD0[i]*PetscPowScalar((1/Vm0),ld_betap[k])*Vr*PetscPowScalar(Vm,(ld_betap[k]-2));
659: dPD_dVi += ld_alphap[k]*ld_betap[k]*PD0[i]*PetscPowScalar((1/Vm0),ld_betap[k])*Vi*PetscPowScalar(Vm,(ld_betap[k]-2));
660: }
661: for (k=0; k < ld_nsegsq[i]; k++) {
662: QD += ld_alphaq[k]*QD0[i]*PetscPowScalar((Vm/Vm0),ld_betaq[k]);
663: dQD_dVr += ld_alphaq[k]*ld_betaq[k]*QD0[i]*PetscPowScalar((1/Vm0),ld_betaq[k])*Vr*PetscPowScalar(Vm,(ld_betaq[k]-2));
664: dQD_dVi += ld_alphaq[k]*ld_betaq[k]*QD0[i]*PetscPowScalar((1/Vm0),ld_betaq[k])*Vi*PetscPowScalar(Vm,(ld_betaq[k]-2));
665: }
667: /* IDr = (PD*Vr + QD*Vi)/Vm2; */
668: /* IDi = (-QD*Vr + PD*Vi)/Vm2; */
670: dIDr_dVr = (dPD_dVr*Vr + dQD_dVr*Vi + PD)/Vm2 - ((PD*Vr + QD*Vi)*2*Vr)/Vm4;
671: dIDr_dVi = (dPD_dVi*Vr + dQD_dVi*Vi + QD)/Vm2 - ((PD*Vr + QD*Vi)*2*Vi)/Vm4;
673: dIDi_dVr = (-dQD_dVr*Vr + dPD_dVr*Vi - QD)/Vm2 - ((-QD*Vr + PD*Vi)*2*Vr)/Vm4;
674: dIDi_dVi = (-dQD_dVi*Vr + dPD_dVi*Vi + PD)/Vm2 - ((-QD*Vr + PD*Vi)*2*Vi)/Vm4;
677: /* fnet[2*lbus[i]] += IDi; */
678: row[0] = net_start + 2*lbus[i];
679: col[0] = net_start + 2*lbus[i]; col[1] = net_start + 2*lbus[i]+1;
680: val[0] = dIDi_dVr; val[1] = dIDi_dVi;
681: MatSetValues(J,1,row,2,col,val,ADD_VALUES);
682: /* fnet[2*lbus[i]+1] += IDr; */
683: row[0] = net_start + 2*lbus[i]+1;
684: col[0] = net_start + 2*lbus[i]; col[1] = net_start + 2*lbus[i]+1;
685: val[0] = dIDr_dVr; val[1] = dIDr_dVi;
686: MatSetValues(J,1,row,2,col,val,ADD_VALUES);
687: }
688: VecRestoreArray(user->V0,&v0);
690: VecRestoreArray(Xgen,&xgen);
691: VecRestoreArray(Xnet,&xnet);
693: DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);
695: MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
696: MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
697: return(0);
698: }
700: /*
701: J = [I, 0
702: dg_dx, dg_dy]
703: */
706: PetscErrorCode AlgJacobian(SNES snes,Vec X,Mat A,Mat B,void *ctx)
707: {
709: Userctx *user=(Userctx*)ctx;
712: ResidualJacobian(snes,X,A,B,ctx);
713: MatSetOption(A,MAT_KEEP_NONZERO_PATTERN,PETSC_TRUE);
714: MatZeroRowsIS(A,user->is_diff,1.0,NULL,NULL);
715: return(0);
716: }
718: /*
719: J = [a*I-df_dx, -df_dy
720: dg_dx, dg_dy]
721: */
725: PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,Userctx *user)
726: {
728: SNES snes;
729: PetscScalar atmp = (PetscScalar) a;
730: PetscInt i,row;
733: user->t = t;
735: TSGetSNES(ts,&snes);
736: ResidualJacobian(snes,X,A,B,user);
737: for (i=0;i < ngen;i++) {
738: row = 9*i;
739: MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
740: row = 9*i+1;
741: MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
742: row = 9*i+2;
743: MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
744: row = 9*i+3;
745: MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
746: row = 9*i+6;
747: MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
748: row = 9*i+7;
749: MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
750: row = 9*i+8;
751: MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
752: }
753: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
754: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
755: return(0);
756: }
760: int main(int argc,char **argv)
761: {
762: TS ts;
763: SNES snes_alg;
765: PetscMPIInt size;
766: Userctx user;
767: PetscViewer Xview,Ybusview;
768: Vec X;
769: Mat J;
770: PetscInt i;
771: /* sensitivity context */
772: PetscScalar *y_ptr;
773: Vec lambda[1];
774: PetscInt steps, total_steps = 0;
775: PetscInt *idx2;
776: Vec Xdot;
777: Vec F_alg;
778: PetscInt row_loc,col_loc;
779: PetscScalar val;
780:
781: PetscInitialize(&argc,&argv,"petscoptions",help);
782: MPI_Comm_size(PETSC_COMM_WORLD,&size);
783: if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs");
785: user.neqs_gen = 9*ngen; /* # eqs. for generator subsystem */
786: user.neqs_net = 2*nbus; /* # eqs. for network subsystem */
787: user.neqs_pgrid = user.neqs_gen + user.neqs_net;
789: /* Create indices for differential and algebraic equations */
790: PetscMalloc1(7*ngen,&idx2);
791: for (i=0; i<ngen; i++) {
792: 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;
793: idx2[7*i+4] = 9*i+6; idx2[7*i+5] = 9*i+7; idx2[7*i+6] = 9*i+8;
794: }
795: ISCreateGeneral(PETSC_COMM_WORLD,7*ngen,idx2,PETSC_COPY_VALUES,&user.is_diff);
796: ISComplement(user.is_diff,0,user.neqs_pgrid,&user.is_alg);
797: PetscFree(idx2);
799: /* Read initial voltage vector and Ybus */
800: PetscViewerBinaryOpen(PETSC_COMM_WORLD,"X.bin",FILE_MODE_READ,&Xview);
801: PetscViewerBinaryOpen(PETSC_COMM_WORLD,"Ybus.bin",FILE_MODE_READ,&Ybusview);
803: VecCreate(PETSC_COMM_WORLD,&user.V0);
804: VecSetSizes(user.V0,PETSC_DECIDE,user.neqs_net);
805: VecLoad(user.V0,Xview);
807: MatCreate(PETSC_COMM_WORLD,&user.Ybus);
808: MatSetSizes(user.Ybus,PETSC_DECIDE,PETSC_DECIDE,user.neqs_net,user.neqs_net);
809: MatSetType(user.Ybus,MATBAIJ);
810: /* MatSetBlockSize(user.Ybus,2); */
811: MatLoad(user.Ybus,Ybusview);
813: /* Set run time options */
814: PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Transient stability fault options","");
815: {
816: user.tfaulton = 1.0;
817: user.tfaultoff = 1.2;
818: user.Rfault = 0.0001;
819: user.faultbus = 8;
820: PetscOptionsReal("-tfaulton","","",user.tfaulton,&user.tfaulton,NULL);
821: PetscOptionsReal("-tfaultoff","","",user.tfaultoff,&user.tfaultoff,NULL);
822: PetscOptionsInt("-faultbus","","",user.faultbus,&user.faultbus,NULL);
823: user.t0 = 0.0;
824: user.tmax = 5.0;
825: PetscOptionsReal("-t0","","",user.t0,&user.t0,NULL);
826: PetscOptionsReal("-tmax","","",user.tmax,&user.tmax,NULL);
827: }
828: PetscOptionsEnd();
830: PetscViewerDestroy(&Xview);
831: PetscViewerDestroy(&Ybusview);
833: /* Create DMs for generator and network subsystems */
834: DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,user.neqs_gen,1,1,NULL,&user.dmgen);
835: DMSetOptionsPrefix(user.dmgen,"dmgen_");
836: DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,user.neqs_net,1,1,NULL,&user.dmnet);
837: DMSetOptionsPrefix(user.dmnet,"dmnet_");
838: /* Create a composite DM packer and add the two DMs */
839: DMCompositeCreate(PETSC_COMM_WORLD,&user.dmpgrid);
840: DMSetOptionsPrefix(user.dmpgrid,"pgrid_");
841: DMCompositeAddDM(user.dmpgrid,user.dmgen);
842: DMCompositeAddDM(user.dmpgrid,user.dmnet);
844: DMCreateGlobalVector(user.dmpgrid,&X);
846: MatCreate(PETSC_COMM_WORLD,&J);
847: MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,user.neqs_pgrid,user.neqs_pgrid);
848: MatSetFromOptions(J);
849: PreallocateJacobian(J,&user);
852: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
853: Create timestepping solver context
854: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
855: TSCreate(PETSC_COMM_WORLD,&ts);
856: TSSetProblemType(ts,TS_NONLINEAR);
857: TSSetEquationType(ts,TS_EQ_DAE_IMPLICIT_INDEX1);
858: TSARKIMEXSetFullyImplicit(ts,PETSC_TRUE);
859: TSSetIFunction(ts,NULL,(TSIFunction) IFunction,&user);
860: TSSetIJacobian(ts,J,J,(TSIJacobian)IJacobian,&user);
861: TSSetApplicationContext(ts,&user);
863: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
864: Set initial conditions
865: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
866: SetInitialGuess(X,&user);
867: /* Just to set up the Jacobian structure */
868: VecDuplicate(X,&Xdot);
869: IJacobian(ts,0.0,X,Xdot,0.0,J,J,&user);
870: VecDestroy(&Xdot);
872: /*
873: Save trajectory of solution so that TSAdjointSolve() may be used
874: */
875: TSSetSaveTrajectory(ts);
877: TSSetDuration(ts,1000,user.tfaulton);
878: TSSetInitialTimeStep(ts,0.0,0.01);
879: TSSetFromOptions(ts);
881: user.alg_flg = PETSC_FALSE;
882: /* Prefault period */
883: TSSolve(ts,X);
884: TSGetTimeStepNumber(ts,&steps);
885: total_steps += steps;
887: /* Create the nonlinear solver for solving the algebraic system */
888: /* Note that although the algebraic system needs to be solved only for
889: Idq and V, we reuse the entire system including xgen. The xgen
890: variables are held constant by setting their residuals to 0 and
891: putting a 1 on the Jacobian diagonal for xgen rows
892: */
893: VecDuplicate(X,&F_alg);
894: SNESCreate(PETSC_COMM_WORLD,&snes_alg);
895: SNESSetFunction(snes_alg,F_alg,AlgFunction,&user);
896: MatZeroEntries(J);
897: SNESSetJacobian(snes_alg,J,J,AlgJacobian,&user);
898: SNESSetOptionsPrefix(snes_alg,"alg_");
899: SNESSetFromOptions(snes_alg);
901: /* Apply disturbance - resistive fault at user.faultbus */
902: /* This is done by adding shunt conductance to the diagonal location
903: in the Ybus matrix */
904: row_loc = 2*user.faultbus; col_loc = 2*user.faultbus+1; /* Location for G */
905: val = 1/user.Rfault;
906: MatSetValues(user.Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);
907: row_loc = 2*user.faultbus+1; col_loc = 2*user.faultbus; /* Location for G */
908: val = 1/user.Rfault;
909: MatSetValues(user.Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);
911: MatAssemblyBegin(user.Ybus,MAT_FINAL_ASSEMBLY);
912: MatAssemblyEnd(user.Ybus,MAT_FINAL_ASSEMBLY);
914: user.alg_flg = PETSC_TRUE;
915: /* Solve the algebraic equations */
916: SNESSolve(snes_alg,NULL,X);
919: /* Disturbance period */
920: TSSetDuration(ts,1000,user.tfaultoff);
921: TSSetInitialTimeStep(ts,user.tfaulton,.01);
923: user.alg_flg = PETSC_FALSE;
925: TSSolve(ts,X);
926: TSGetTimeStepNumber(ts,&steps);
927: total_steps += steps;
928: /* Remove the fault */
929: row_loc = 2*user.faultbus; col_loc = 2*user.faultbus+1;
930: val = -1/user.Rfault;
931: MatSetValues(user.Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);
932: row_loc = 2*user.faultbus+1; col_loc = 2*user.faultbus;
933: val = -1/user.Rfault;
934: MatSetValues(user.Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);
936: MatAssemblyBegin(user.Ybus,MAT_FINAL_ASSEMBLY);
937: MatAssemblyEnd(user.Ybus,MAT_FINAL_ASSEMBLY);
939: MatZeroEntries(J);
941: user.alg_flg = PETSC_TRUE;
943: /* Solve the algebraic equations */
944: SNESSolve(snes_alg,NULL,X);
946: /* Post-disturbance period */
947: TSSetDuration(ts,1000,user.tmax);
948: TSSetInitialTimeStep(ts,user.tfaultoff,.01);
950: user.alg_flg = PETSC_TRUE;
952: TSSolve(ts,X);
953: TSGetTimeStepNumber(ts,&steps);
954: total_steps += steps;
956: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
957: Adjoint model starts here
958: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
959: TSSetPostStep(ts,NULL);
960: MatCreateVecs(J,&lambda[0],NULL);
961: /* Set initial conditions for the adjoint integration */
962: VecZeroEntries(lambda[0]);
963: VecGetArray(lambda[0],&y_ptr);
964: y_ptr[0] = 1.0;
965: VecRestoreArray(lambda[0],&y_ptr);
966: TSSetCostGradients(ts,1,lambda,NULL);
968: TSAdjointSolve(ts);
970: PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt initial conditions: \n");
971: VecView(lambda[0],PETSC_VIEWER_STDOUT_WORLD);
972: VecDestroy(&lambda[0]);
974: SNESDestroy(&snes_alg);
975: VecDestroy(&F_alg);
976: MatDestroy(&J);
977: MatDestroy(&user.Ybus);
978: VecDestroy(&X);
979: VecDestroy(&user.V0);
980: DMDestroy(&user.dmgen);
981: DMDestroy(&user.dmnet);
982: DMDestroy(&user.dmpgrid);
983: ISDestroy(&user.is_diff);
984: ISDestroy(&user.is_alg);
985: TSDestroy(&ts);
986: PetscFinalize();
987: return(0);
988: }