xref: /petsc/src/ts/tutorials/power_grid/stability_9bus/ex9busopt.c (revision 6a98f8dc3f2c9149905a87dc2e9d0fedaf64e09a)
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";
7 
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 */
14 
15 #include <petsctao.h>
16 #include <petscts.h>
17 #include <petscdm.h>
18 #include <petscdmda.h>
19 #include <petscdmcomposite.h>
20 #include <petsctime.h>
21 
22 PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal*,Vec,void*);
23 
24 #define freq 60
25 #define w_s (2*PETSC_PI*freq)
26 
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 */
33 
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};*/
37 
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 */
50 
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) */
61 
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
75 
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};
86 
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   Mat         DRDU,DRDP;
108 } Userctx;
109 
110 
111 /* Converts from machine frame (dq) to network (phase a real,imag) reference frame */
112 PetscErrorCode dq2ri(PetscScalar Fd,PetscScalar Fq,PetscScalar delta,PetscScalar *Fr, PetscScalar *Fi)
113 {
114   PetscFunctionBegin;
115   *Fr =  Fd*PetscSinScalar(delta) + Fq*PetscCosScalar(delta);
116   *Fi = -Fd*PetscCosScalar(delta) + Fq*PetscSinScalar(delta);
117   PetscFunctionReturn(0);
118 }
119 
120 /* Converts from network frame ([phase a real,imag) to machine (dq) reference frame */
121 PetscErrorCode ri2dq(PetscScalar Fr,PetscScalar Fi,PetscScalar delta,PetscScalar *Fd, PetscScalar *Fq)
122 {
123   PetscFunctionBegin;
124   *Fd =  Fr*PetscSinScalar(delta) - Fi*PetscCosScalar(delta);
125   *Fq =  Fr*PetscCosScalar(delta) + Fi*PetscSinScalar(delta);
126   PetscFunctionReturn(0);
127 }
128 
129 /* Saves the solution at each time to a matrix */
130 PetscErrorCode SaveSolution(TS ts)
131 {
132   PetscErrorCode    ierr;
133   Userctx           *user;
134   Vec               X;
135   PetscScalar       *mat;
136   const PetscScalar *x;
137   PetscInt          idx;
138   PetscReal         t;
139 
140   PetscFunctionBegin;
141   ierr     = TSGetApplicationContext(ts,&user);CHKERRQ(ierr);
142   ierr     = TSGetTime(ts,&t);CHKERRQ(ierr);
143   ierr     = TSGetSolution(ts,&X);CHKERRQ(ierr);
144   idx      = user->stepnum*(user->neqs_pgrid+1);
145   ierr     = MatDenseGetArray(user->Sol,&mat);CHKERRQ(ierr);
146   ierr     = VecGetArrayRead(X,&x);CHKERRQ(ierr);
147   mat[idx] = t;
148   ierr     = PetscArraycpy(mat+idx+1,x,user->neqs_pgrid);CHKERRQ(ierr);
149   ierr     = MatDenseRestoreArray(user->Sol,&mat);CHKERRQ(ierr);
150   ierr     = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
151   user->stepnum++;
152   PetscFunctionReturn(0);
153 }
154 
155 PetscErrorCode SetInitialGuess(Vec X,Userctx *user)
156 {
157   PetscErrorCode ierr;
158   Vec            Xgen,Xnet;
159   PetscScalar    *xgen,*xnet;
160   PetscInt       i,idx=0;
161   PetscScalar    Vr,Vi,IGr,IGi,Vm,Vm2;
162   PetscScalar    Eqp,Edp,delta;
163   PetscScalar    Efd,RF,VR; /* Exciter variables */
164   PetscScalar    Id,Iq;  /* Generator dq axis currents */
165   PetscScalar    theta,Vd,Vq,SE;
166 
167   PetscFunctionBegin;
168   M[0] = 2*H[0]/w_s; M[1] = 2*H[1]/w_s; M[2] = 2*H[2]/w_s;
169   D[0] = 0.1*M[0]; D[1] = 0.1*M[1]; D[2] = 0.1*M[2];
170 
171   ierr = DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);CHKERRQ(ierr);
172 
173   /* Network subsystem initialization */
174   ierr = VecCopy(user->V0,Xnet);CHKERRQ(ierr);
175 
176   /* Generator subsystem initialization */
177   ierr = VecGetArray(Xgen,&xgen);CHKERRQ(ierr);
178   ierr = VecGetArray(Xnet,&xnet);CHKERRQ(ierr);
179 
180   for (i=0; i < ngen; i++) {
181     Vr  = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */
182     Vi  = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */
183     Vm  = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm;
184     IGr = (Vr*PG[i] + Vi*QG[i])/Vm2;
185     IGi = (Vi*PG[i] - Vr*QG[i])/Vm2;
186 
187     delta = PetscAtan2Real(Vi+Xq[i]*IGr,Vr-Xq[i]*IGi); /* Machine angle */
188 
189     theta = PETSC_PI/2.0 - delta;
190 
191     Id = IGr*PetscCosScalar(theta) - IGi*PetscSinScalar(theta); /* d-axis stator current */
192     Iq = IGr*PetscSinScalar(theta) + IGi*PetscCosScalar(theta); /* q-axis stator current */
193 
194     Vd = Vr*PetscCosScalar(theta) - Vi*PetscSinScalar(theta);
195     Vq = Vr*PetscSinScalar(theta) + Vi*PetscCosScalar(theta);
196 
197     Edp = Vd + Rs[i]*Id - Xqp[i]*Iq; /* d-axis transient EMF */
198     Eqp = Vq + Rs[i]*Iq + Xdp[i]*Id; /* q-axis transient EMF */
199 
200     TM[i] = PG[i];
201 
202     /* The generator variables are ordered as [Eqp,Edp,delta,w,Id,Iq] */
203     xgen[idx]   = Eqp;
204     xgen[idx+1] = Edp;
205     xgen[idx+2] = delta;
206     xgen[idx+3] = w_s;
207 
208     idx = idx + 4;
209 
210     xgen[idx]   = Id;
211     xgen[idx+1] = Iq;
212 
213     idx = idx + 2;
214 
215     /* Exciter */
216     Efd = Eqp + (Xd[i] - Xdp[i])*Id;
217     SE  = k1[i]*PetscExpScalar(k2[i]*Efd);
218     VR  =  KE[i]*Efd + SE;
219     RF  =  KF[i]*Efd/TF[i];
220 
221     xgen[idx]   = Efd;
222     xgen[idx+1] = RF;
223     xgen[idx+2] = VR;
224 
225     Vref[i] = Vm + (VR/KA[i]);
226 
227     idx = idx + 3;
228   }
229 
230   ierr = VecRestoreArray(Xgen,&xgen);CHKERRQ(ierr);
231   ierr = VecRestoreArray(Xnet,&xnet);CHKERRQ(ierr);
232 
233   /* ierr = VecView(Xgen,0);CHKERRQ(ierr); */
234   ierr = DMCompositeGather(user->dmpgrid,INSERT_VALUES,X,Xgen,Xnet);CHKERRQ(ierr);
235   ierr = DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);CHKERRQ(ierr);
236   PetscFunctionReturn(0);
237 }
238 
239 PetscErrorCode InitialGuess(Vec X,Userctx *user, const PetscScalar PGv[])
240 {
241   PetscErrorCode ierr;
242   Vec            Xgen,Xnet;
243   PetscScalar    *xgen,*xnet;
244   PetscInt       i,idx=0;
245   PetscScalar    Vr,Vi,IGr,IGi,Vm,Vm2;
246   PetscScalar    Eqp,Edp,delta;
247   PetscScalar    Efd,RF,VR; /* Exciter variables */
248   PetscScalar    Id,Iq;  /* Generator dq axis currents */
249   PetscScalar    theta,Vd,Vq,SE;
250 
251   PetscFunctionBegin;
252   M[0] = 2*H[0]/w_s; M[1] = 2*H[1]/w_s; M[2] = 2*H[2]/w_s;
253   D[0] = 0.1*M[0]; D[1] = 0.1*M[1]; D[2] = 0.1*M[2];
254 
255   ierr = DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);CHKERRQ(ierr);
256 
257   /* Network subsystem initialization */
258   ierr = VecCopy(user->V0,Xnet);CHKERRQ(ierr);
259 
260   /* Generator subsystem initialization */
261   ierr = VecGetArray(Xgen,&xgen);CHKERRQ(ierr);
262   ierr = VecGetArray(Xnet,&xnet);CHKERRQ(ierr);
263 
264   for (i=0; i < ngen; i++) {
265     Vr  = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */
266     Vi  = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */
267     Vm  = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm;
268     IGr = (Vr*PGv[i] + Vi*QG[i])/Vm2;
269     IGi = (Vi*PGv[i] - Vr*QG[i])/Vm2;
270 
271     delta = PetscAtan2Real(Vi+Xq[i]*IGr,Vr-Xq[i]*IGi); /* Machine angle */
272 
273     theta = PETSC_PI/2.0 - delta;
274 
275     Id = IGr*PetscCosScalar(theta) - IGi*PetscSinScalar(theta); /* d-axis stator current */
276     Iq = IGr*PetscSinScalar(theta) + IGi*PetscCosScalar(theta); /* q-axis stator current */
277 
278     Vd = Vr*PetscCosScalar(theta) - Vi*PetscSinScalar(theta);
279     Vq = Vr*PetscSinScalar(theta) + Vi*PetscCosScalar(theta);
280 
281     Edp = Vd + Rs[i]*Id - Xqp[i]*Iq; /* d-axis transient EMF */
282     Eqp = Vq + Rs[i]*Iq + Xdp[i]*Id; /* q-axis transient EMF */
283 
284     /* The generator variables are ordered as [Eqp,Edp,delta,w,Id,Iq] */
285     xgen[idx]   = Eqp;
286     xgen[idx+1] = Edp;
287     xgen[idx+2] = delta;
288     xgen[idx+3] = w_s;
289 
290     idx = idx + 4;
291 
292     xgen[idx]   = Id;
293     xgen[idx+1] = Iq;
294 
295     idx = idx + 2;
296 
297     /* Exciter */
298     Efd = Eqp + (Xd[i] - Xdp[i])*Id;
299     SE  = k1[i]*PetscExpScalar(k2[i]*Efd);
300     VR  =  KE[i]*Efd + SE;
301     RF  =  KF[i]*Efd/TF[i];
302 
303     xgen[idx]   = Efd;
304     xgen[idx+1] = RF;
305     xgen[idx+2] = VR;
306 
307     idx = idx + 3;
308   }
309 
310   ierr = VecRestoreArray(Xgen,&xgen);CHKERRQ(ierr);
311   ierr = VecRestoreArray(Xnet,&xnet);CHKERRQ(ierr);
312 
313   /* ierr = VecView(Xgen,0);CHKERRQ(ierr); */
314   ierr = DMCompositeGather(user->dmpgrid,INSERT_VALUES,X,Xgen,Xnet);CHKERRQ(ierr);
315   ierr = DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);CHKERRQ(ierr);
316   PetscFunctionReturn(0);
317 }
318 
319 PetscErrorCode DICDPFiniteDifference(Vec X,Vec *DICDP, Userctx *user)
320 {
321   Vec            Y;
322   PetscScalar    PGv[3],eps;
323   PetscErrorCode ierr;
324   PetscInt       i,j;
325 
326   eps = 1.e-7;
327   ierr = VecDuplicate(X,&Y);CHKERRQ(ierr);
328 
329   for (i=0;i<ngen;i++) {
330     for (j=0;j<3;j++) PGv[j] = PG[j];
331     PGv[i] = PG[i]+eps;
332     ierr = InitialGuess(Y,user,PGv);CHKERRQ(ierr);
333     ierr = InitialGuess(X,user,PG);CHKERRQ(ierr);
334 
335     ierr = VecAXPY(Y,-1.0,X);CHKERRQ(ierr);
336     ierr = VecScale(Y,1./eps);CHKERRQ(ierr);
337     ierr = VecCopy(Y,DICDP[i]);CHKERRQ(ierr);
338   }
339   ierr = VecDestroy(&Y);CHKERRQ(ierr);
340   PetscFunctionReturn(0);
341 }
342 
343 
344 /* Computes F = [-f(x,y);g(x,y)] */
345 PetscErrorCode ResidualFunction(SNES snes,Vec X, Vec F, Userctx *user)
346 {
347   PetscErrorCode ierr;
348   Vec            Xgen,Xnet,Fgen,Fnet;
349   PetscScalar    *xgen,*xnet,*fgen,*fnet;
350   PetscInt       i,idx=0;
351   PetscScalar    Vr,Vi,Vm,Vm2;
352   PetscScalar    Eqp,Edp,delta,w; /* Generator variables */
353   PetscScalar    Efd,RF,VR; /* Exciter variables */
354   PetscScalar    Id,Iq;  /* Generator dq axis currents */
355   PetscScalar    Vd,Vq,SE;
356   PetscScalar    IGr,IGi,IDr,IDi;
357   PetscScalar    Zdq_inv[4],det;
358   PetscScalar    PD,QD,Vm0,*v0;
359   PetscInt       k;
360 
361   PetscFunctionBegin;
362   ierr = VecZeroEntries(F);CHKERRQ(ierr);
363   ierr = DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);CHKERRQ(ierr);
364   ierr = DMCompositeGetLocalVectors(user->dmpgrid,&Fgen,&Fnet);CHKERRQ(ierr);
365   ierr = DMCompositeScatter(user->dmpgrid,X,Xgen,Xnet);CHKERRQ(ierr);
366   ierr = DMCompositeScatter(user->dmpgrid,F,Fgen,Fnet);CHKERRQ(ierr);
367 
368   /* Network current balance residual IG + Y*V + IL = 0. Only YV is added here.
369      The generator current injection, IG, and load current injection, ID are added later
370   */
371   /* Note that the values in Ybus are stored assuming the imaginary current balance
372      equation is ordered first followed by real current balance equation for each bus.
373      Thus imaginary current contribution goes in location 2*i, and
374      real current contribution in 2*i+1
375   */
376   ierr = MatMult(user->Ybus,Xnet,Fnet);CHKERRQ(ierr);
377 
378   ierr = VecGetArray(Xgen,&xgen);CHKERRQ(ierr);
379   ierr = VecGetArray(Xnet,&xnet);CHKERRQ(ierr);
380   ierr = VecGetArray(Fgen,&fgen);CHKERRQ(ierr);
381   ierr = VecGetArray(Fnet,&fnet);CHKERRQ(ierr);
382 
383   /* Generator subsystem */
384   for (i=0; i < ngen; i++) {
385     Eqp   = xgen[idx];
386     Edp   = xgen[idx+1];
387     delta = xgen[idx+2];
388     w     = xgen[idx+3];
389     Id    = xgen[idx+4];
390     Iq    = xgen[idx+5];
391     Efd   = xgen[idx+6];
392     RF    = xgen[idx+7];
393     VR    = xgen[idx+8];
394 
395     /* Generator differential equations */
396     fgen[idx]   = (Eqp + (Xd[i] - Xdp[i])*Id - Efd)/Td0p[i];
397     fgen[idx+1] = (Edp - (Xq[i] - Xqp[i])*Iq)/Tq0p[i];
398     fgen[idx+2] = -w + w_s;
399     fgen[idx+3] = (-TM[i] + Edp*Id + Eqp*Iq + (Xqp[i] - Xdp[i])*Id*Iq + D[i]*(w - w_s))/M[i];
400 
401     Vr = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */
402     Vi = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */
403 
404     ierr = ri2dq(Vr,Vi,delta,&Vd,&Vq);CHKERRQ(ierr);
405     /* Algebraic equations for stator currents */
406     det = Rs[i]*Rs[i] + Xdp[i]*Xqp[i];
407 
408     Zdq_inv[0] = Rs[i]/det;
409     Zdq_inv[1] = Xqp[i]/det;
410     Zdq_inv[2] = -Xdp[i]/det;
411     Zdq_inv[3] = Rs[i]/det;
412 
413     fgen[idx+4] = Zdq_inv[0]*(-Edp + Vd) + Zdq_inv[1]*(-Eqp + Vq) + Id;
414     fgen[idx+5] = Zdq_inv[2]*(-Edp + Vd) + Zdq_inv[3]*(-Eqp + Vq) + Iq;
415 
416     /* Add generator current injection to network */
417     ierr = dq2ri(Id,Iq,delta,&IGr,&IGi);CHKERRQ(ierr);
418 
419     fnet[2*gbus[i]]   -= IGi;
420     fnet[2*gbus[i]+1] -= IGr;
421 
422     Vm = PetscSqrtScalar(Vd*Vd + Vq*Vq);
423 
424     SE = k1[i]*PetscExpScalar(k2[i]*Efd);
425 
426     /* Exciter differential equations */
427     fgen[idx+6] = (KE[i]*Efd + SE - VR)/TE[i];
428     fgen[idx+7] = (RF - KF[i]*Efd/TF[i])/TF[i];
429     fgen[idx+8] = (VR - KA[i]*RF + KA[i]*KF[i]*Efd/TF[i] - KA[i]*(Vref[i] - Vm))/TA[i];
430 
431     idx = idx + 9;
432   }
433 
434   ierr = VecGetArray(user->V0,&v0);CHKERRQ(ierr);
435   for (i=0; i < nload; i++) {
436     Vr  = xnet[2*lbus[i]]; /* Real part of load bus voltage */
437     Vi  = xnet[2*lbus[i]+1]; /* Imaginary part of the load bus voltage */
438     Vm  = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm;
439     Vm0 = PetscSqrtScalar(v0[2*lbus[i]]*v0[2*lbus[i]] + v0[2*lbus[i]+1]*v0[2*lbus[i]+1]);
440     PD  = QD = 0.0;
441     for (k=0; k < ld_nsegsp[i]; k++) PD += ld_alphap[k]*PD0[i]*PetscPowScalar((Vm/Vm0),ld_betap[k]);
442     for (k=0; k < ld_nsegsq[i]; k++) QD += ld_alphaq[k]*QD0[i]*PetscPowScalar((Vm/Vm0),ld_betaq[k]);
443 
444     /* Load currents */
445     IDr = (PD*Vr + QD*Vi)/Vm2;
446     IDi = (-QD*Vr + PD*Vi)/Vm2;
447 
448     fnet[2*lbus[i]]   += IDi;
449     fnet[2*lbus[i]+1] += IDr;
450   }
451   ierr = VecRestoreArray(user->V0,&v0);CHKERRQ(ierr);
452 
453   ierr = VecRestoreArray(Xgen,&xgen);CHKERRQ(ierr);
454   ierr = VecRestoreArray(Xnet,&xnet);CHKERRQ(ierr);
455   ierr = VecRestoreArray(Fgen,&fgen);CHKERRQ(ierr);
456   ierr = VecRestoreArray(Fnet,&fnet);CHKERRQ(ierr);
457 
458   ierr = DMCompositeGather(user->dmpgrid,INSERT_VALUES,F,Fgen,Fnet);CHKERRQ(ierr);
459   ierr = DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);CHKERRQ(ierr);
460   ierr = DMCompositeRestoreLocalVectors(user->dmpgrid,&Fgen,&Fnet);CHKERRQ(ierr);
461   PetscFunctionReturn(0);
462 }
463 
464 /* \dot{x} - f(x,y)
465      g(x,y) = 0
466  */
467 PetscErrorCode IFunction(TS ts,PetscReal t, Vec X, Vec Xdot, Vec F, Userctx *user)
468 {
469   PetscErrorCode    ierr;
470   SNES              snes;
471   PetscScalar       *f;
472   const PetscScalar *xdot;
473   PetscInt          i;
474 
475   PetscFunctionBegin;
476   user->t = t;
477 
478   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
479   ierr = ResidualFunction(snes,X,F,user);CHKERRQ(ierr);
480   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
481   ierr = VecGetArrayRead(Xdot,&xdot);CHKERRQ(ierr);
482   for (i=0;i < ngen;i++) {
483     f[9*i]   += xdot[9*i];
484     f[9*i+1] += xdot[9*i+1];
485     f[9*i+2] += xdot[9*i+2];
486     f[9*i+3] += xdot[9*i+3];
487     f[9*i+6] += xdot[9*i+6];
488     f[9*i+7] += xdot[9*i+7];
489     f[9*i+8] += xdot[9*i+8];
490   }
491   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
492   ierr = VecRestoreArrayRead(Xdot,&xdot);CHKERRQ(ierr);
493   PetscFunctionReturn(0);
494 }
495 
496 /* This function is used for solving the algebraic system only during fault on and
497    off times. It computes the entire F and then zeros out the part corresponding to
498    differential equations
499  F = [0;g(y)];
500 */
501 PetscErrorCode AlgFunction(SNES snes, Vec X, Vec F, void *ctx)
502 {
503   PetscErrorCode ierr;
504   Userctx        *user=(Userctx*)ctx;
505   PetscScalar    *f;
506   PetscInt       i;
507 
508   PetscFunctionBegin;
509   ierr = ResidualFunction(snes,X,F,user);CHKERRQ(ierr);
510   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
511   for (i=0; i < ngen; i++) {
512     f[9*i]   = 0;
513     f[9*i+1] = 0;
514     f[9*i+2] = 0;
515     f[9*i+3] = 0;
516     f[9*i+6] = 0;
517     f[9*i+7] = 0;
518     f[9*i+8] = 0;
519   }
520   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
521   PetscFunctionReturn(0);
522 }
523 
524 PetscErrorCode PreallocateJacobian(Mat J, Userctx *user)
525 {
526   PetscErrorCode ierr;
527   PetscInt       *d_nnz;
528   PetscInt       i,idx=0,start=0;
529   PetscInt       ncols;
530 
531   PetscFunctionBegin;
532   ierr = PetscMalloc1(user->neqs_pgrid,&d_nnz);CHKERRQ(ierr);
533   for (i=0; i<user->neqs_pgrid; i++) d_nnz[i] = 0;
534   /* Generator subsystem */
535   for (i=0; i < ngen; i++) {
536 
537     d_nnz[idx]   += 3;
538     d_nnz[idx+1] += 2;
539     d_nnz[idx+2] += 2;
540     d_nnz[idx+3] += 5;
541     d_nnz[idx+4] += 6;
542     d_nnz[idx+5] += 6;
543 
544     d_nnz[user->neqs_gen+2*gbus[i]]   += 3;
545     d_nnz[user->neqs_gen+2*gbus[i]+1] += 3;
546 
547     d_nnz[idx+6] += 2;
548     d_nnz[idx+7] += 2;
549     d_nnz[idx+8] += 5;
550 
551     idx = idx + 9;
552   }
553 
554   start = user->neqs_gen;
555   for (i=0; i < nbus; i++) {
556     ierr = MatGetRow(user->Ybus,2*i,&ncols,NULL,NULL);CHKERRQ(ierr);
557     d_nnz[start+2*i]   += ncols;
558     d_nnz[start+2*i+1] += ncols;
559     ierr = MatRestoreRow(user->Ybus,2*i,&ncols,NULL,NULL);CHKERRQ(ierr);
560   }
561 
562   ierr = MatSeqAIJSetPreallocation(J,0,d_nnz);CHKERRQ(ierr);
563   ierr = PetscFree(d_nnz);CHKERRQ(ierr);
564   PetscFunctionReturn(0);
565 }
566 
567 /*
568    J = [-df_dx, -df_dy
569         dg_dx, dg_dy]
570 */
571 PetscErrorCode ResidualJacobian(SNES snes,Vec X,Mat J,Mat B,void *ctx)
572 {
573   PetscErrorCode    ierr;
574   Userctx           *user=(Userctx*)ctx;
575   Vec               Xgen,Xnet;
576   PetscScalar       *xgen,*xnet;
577   PetscInt          i,idx=0;
578   PetscScalar       Vr,Vi,Vm,Vm2;
579   PetscScalar       Eqp,Edp,delta; /* Generator variables */
580   PetscScalar       Efd; /* Exciter variables */
581   PetscScalar       Id,Iq;  /* Generator dq axis currents */
582   PetscScalar       Vd,Vq;
583   PetscScalar       val[10];
584   PetscInt          row[2],col[10];
585   PetscInt          net_start=user->neqs_gen;
586   PetscInt          ncols;
587   const PetscInt    *cols;
588   const PetscScalar *yvals;
589   PetscInt          k;
590   PetscScalar       Zdq_inv[4],det;
591   PetscScalar       dVd_dVr,dVd_dVi,dVq_dVr,dVq_dVi,dVd_ddelta,dVq_ddelta;
592   PetscScalar       dIGr_ddelta,dIGi_ddelta,dIGr_dId,dIGr_dIq,dIGi_dId,dIGi_dIq;
593   PetscScalar       dSE_dEfd;
594   PetscScalar       dVm_dVd,dVm_dVq,dVm_dVr,dVm_dVi;
595   PetscScalar       PD,QD,Vm0,*v0,Vm4;
596   PetscScalar       dPD_dVr,dPD_dVi,dQD_dVr,dQD_dVi;
597   PetscScalar       dIDr_dVr,dIDr_dVi,dIDi_dVr,dIDi_dVi;
598 
599   PetscFunctionBegin;
600   ierr  = MatZeroEntries(B);CHKERRQ(ierr);
601   ierr  = DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);CHKERRQ(ierr);
602   ierr  = DMCompositeScatter(user->dmpgrid,X,Xgen,Xnet);CHKERRQ(ierr);
603 
604   ierr = VecGetArray(Xgen,&xgen);CHKERRQ(ierr);
605   ierr = VecGetArray(Xnet,&xnet);CHKERRQ(ierr);
606 
607   /* Generator subsystem */
608   for (i=0; i < ngen; i++) {
609     Eqp   = xgen[idx];
610     Edp   = xgen[idx+1];
611     delta = xgen[idx+2];
612     Id    = xgen[idx+4];
613     Iq    = xgen[idx+5];
614     Efd   = xgen[idx+6];
615 
616     /*    fgen[idx]   = (Eqp + (Xd[i] - Xdp[i])*Id - Efd)/Td0p[i]; */
617     row[0] = idx;
618     col[0] = idx;           col[1] = idx+4;          col[2] = idx+6;
619     val[0] = 1/ Td0p[i]; val[1] = (Xd[i] - Xdp[i])/ Td0p[i]; val[2] = -1/Td0p[i];
620 
621     ierr = MatSetValues(J,1,row,3,col,val,INSERT_VALUES);CHKERRQ(ierr);
622 
623     /*    fgen[idx+1] = (Edp - (Xq[i] - Xqp[i])*Iq)/Tq0p[i]; */
624     row[0] = idx + 1;
625     col[0] = idx + 1;       col[1] = idx+5;
626     val[0] = 1/Tq0p[i]; val[1] = -(Xq[i] - Xqp[i])/Tq0p[i];
627     ierr   = MatSetValues(J,1,row,2,col,val,INSERT_VALUES);CHKERRQ(ierr);
628 
629     /*    fgen[idx+2] = - w + w_s; */
630     row[0] = idx + 2;
631     col[0] = idx + 2; col[1] = idx + 3;
632     val[0] = 0;       val[1] = -1;
633     ierr   = MatSetValues(J,1,row,2,col,val,INSERT_VALUES);CHKERRQ(ierr);
634 
635     /*    fgen[idx+3] = (-TM[i] + Edp*Id + Eqp*Iq + (Xqp[i] - Xdp[i])*Id*Iq + D[i]*(w - w_s))/M[i]; */
636     row[0] = idx + 3;
637     col[0] = idx; col[1] = idx + 1; col[2] = idx + 3;       col[3] = idx + 4;                  col[4] = idx + 5;
638     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];
639     ierr   = MatSetValues(J,1,row,5,col,val,INSERT_VALUES);CHKERRQ(ierr);
640 
641     Vr   = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */
642     Vi   = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */
643     ierr = ri2dq(Vr,Vi,delta,&Vd,&Vq);CHKERRQ(ierr);
644 
645     det = Rs[i]*Rs[i] + Xdp[i]*Xqp[i];
646 
647     Zdq_inv[0] = Rs[i]/det;
648     Zdq_inv[1] = Xqp[i]/det;
649     Zdq_inv[2] = -Xdp[i]/det;
650     Zdq_inv[3] = Rs[i]/det;
651 
652     dVd_dVr    = PetscSinScalar(delta); dVd_dVi = -PetscCosScalar(delta);
653     dVq_dVr    = PetscCosScalar(delta); dVq_dVi = PetscSinScalar(delta);
654     dVd_ddelta = Vr*PetscCosScalar(delta) + Vi*PetscSinScalar(delta);
655     dVq_ddelta = -Vr*PetscSinScalar(delta) + Vi*PetscCosScalar(delta);
656 
657     /*    fgen[idx+4] = Zdq_inv[0]*(-Edp + Vd) + Zdq_inv[1]*(-Eqp + Vq) + Id; */
658     row[0] = idx+4;
659     col[0] = idx;         col[1] = idx+1;        col[2] = idx + 2;
660     val[0] = -Zdq_inv[1]; val[1] = -Zdq_inv[0];  val[2] = Zdq_inv[0]*dVd_ddelta + Zdq_inv[1]*dVq_ddelta;
661     col[3] = idx + 4; col[4] = net_start+2*gbus[i];                     col[5] = net_start + 2*gbus[i]+1;
662     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;
663     ierr   = MatSetValues(J,1,row,6,col,val,INSERT_VALUES);CHKERRQ(ierr);
664 
665     /*  fgen[idx+5] = Zdq_inv[2]*(-Edp + Vd) + Zdq_inv[3]*(-Eqp + Vq) + Iq; */
666     row[0] = idx+5;
667     col[0] = idx;         col[1] = idx+1;        col[2] = idx + 2;
668     val[0] = -Zdq_inv[3]; val[1] = -Zdq_inv[2];  val[2] = Zdq_inv[2]*dVd_ddelta + Zdq_inv[3]*dVq_ddelta;
669     col[3] = idx + 5; col[4] = net_start+2*gbus[i];                     col[5] = net_start + 2*gbus[i]+1;
670     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;
671     ierr   = MatSetValues(J,1,row,6,col,val,INSERT_VALUES);CHKERRQ(ierr);
672 
673     dIGr_ddelta = Id*PetscCosScalar(delta) - Iq*PetscSinScalar(delta);
674     dIGi_ddelta = Id*PetscSinScalar(delta) + Iq*PetscCosScalar(delta);
675     dIGr_dId    = PetscSinScalar(delta);  dIGr_dIq = PetscCosScalar(delta);
676     dIGi_dId    = -PetscCosScalar(delta); dIGi_dIq = PetscSinScalar(delta);
677 
678     /* fnet[2*gbus[i]]   -= IGi; */
679     row[0] = net_start + 2*gbus[i];
680     col[0] = idx+2;        col[1] = idx + 4;   col[2] = idx + 5;
681     val[0] = -dIGi_ddelta; val[1] = -dIGi_dId; val[2] = -dIGi_dIq;
682     ierr = MatSetValues(J,1,row,3,col,val,INSERT_VALUES);CHKERRQ(ierr);
683 
684     /* fnet[2*gbus[i]+1]   -= IGr; */
685     row[0] = net_start + 2*gbus[i]+1;
686     col[0] = idx+2;        col[1] = idx + 4;   col[2] = idx + 5;
687     val[0] = -dIGr_ddelta; val[1] = -dIGr_dId; val[2] = -dIGr_dIq;
688     ierr   = MatSetValues(J,1,row,3,col,val,INSERT_VALUES);CHKERRQ(ierr);
689 
690     Vm = PetscSqrtScalar(Vd*Vd + Vq*Vq);
691 
692     /*    fgen[idx+6] = (KE[i]*Efd + SE - VR)/TE[i]; */
693     /*    SE  = k1[i]*PetscExpScalar(k2[i]*Efd); */
694     dSE_dEfd = k1[i]*k2[i]*PetscExpScalar(k2[i]*Efd);
695 
696     row[0] = idx + 6;
697     col[0] = idx + 6;                     col[1] = idx + 8;
698     val[0] = (KE[i] + dSE_dEfd)/TE[i];  val[1] = -1/TE[i];
699     ierr   = MatSetValues(J,1,row,2,col,val,INSERT_VALUES);CHKERRQ(ierr);
700 
701     /* Exciter differential equations */
702 
703     /*    fgen[idx+7] = (RF - KF[i]*Efd/TF[i])/TF[i]; */
704     row[0] = idx + 7;
705     col[0] = idx + 6;       col[1] = idx + 7;
706     val[0] = (-KF[i]/TF[i])/TF[i];  val[1] = 1/TF[i];
707     ierr   = MatSetValues(J,1,row,2,col,val,INSERT_VALUES);CHKERRQ(ierr);
708 
709     /*    fgen[idx+8] = (VR - KA[i]*RF + KA[i]*KF[i]*Efd/TF[i] - KA[i]*(Vref[i] - Vm))/TA[i]; */
710     /* Vm = (Vd^2 + Vq^2)^0.5; */
711     dVm_dVd    = Vd/Vm; dVm_dVq = Vq/Vm;
712     dVm_dVr    = dVm_dVd*dVd_dVr + dVm_dVq*dVq_dVr;
713     dVm_dVi    = dVm_dVd*dVd_dVi + dVm_dVq*dVq_dVi;
714     row[0]     = idx + 8;
715     col[0]     = idx + 6;           col[1] = idx + 7; col[2] = idx + 8;
716     val[0]     = (KA[i]*KF[i]/TF[i])/TA[i]; val[1] = -KA[i]/TA[i];  val[2] = 1/TA[i];
717     col[3]     = net_start + 2*gbus[i]; col[4] = net_start + 2*gbus[i]+1;
718     val[3]     = KA[i]*dVm_dVr/TA[i];         val[4] = KA[i]*dVm_dVi/TA[i];
719     ierr       = MatSetValues(J,1,row,5,col,val,INSERT_VALUES);CHKERRQ(ierr);
720     idx        = idx + 9;
721   }
722 
723 
724   for (i=0; i<nbus; i++) {
725     ierr   = MatGetRow(user->Ybus,2*i,&ncols,&cols,&yvals);CHKERRQ(ierr);
726     row[0] = net_start + 2*i;
727     for (k=0; k<ncols; k++) {
728       col[k] = net_start + cols[k];
729       val[k] = yvals[k];
730     }
731     ierr = MatSetValues(J,1,row,ncols,col,val,INSERT_VALUES);CHKERRQ(ierr);
732     ierr = MatRestoreRow(user->Ybus,2*i,&ncols,&cols,&yvals);CHKERRQ(ierr);
733 
734     ierr   = MatGetRow(user->Ybus,2*i+1,&ncols,&cols,&yvals);CHKERRQ(ierr);
735     row[0] = net_start + 2*i+1;
736     for (k=0; k<ncols; k++) {
737       col[k] = net_start + cols[k];
738       val[k] = yvals[k];
739     }
740     ierr = MatSetValues(J,1,row,ncols,col,val,INSERT_VALUES);CHKERRQ(ierr);
741     ierr = MatRestoreRow(user->Ybus,2*i+1,&ncols,&cols,&yvals);CHKERRQ(ierr);
742   }
743 
744   ierr = MatAssemblyBegin(J,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr);
745   ierr = MatAssemblyEnd(J,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr);
746 
747 
748   ierr = VecGetArray(user->V0,&v0);CHKERRQ(ierr);
749   for (i=0; i < nload; i++) {
750     Vr      = xnet[2*lbus[i]]; /* Real part of load bus voltage */
751     Vi      = xnet[2*lbus[i]+1]; /* Imaginary part of the load bus voltage */
752     Vm      = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2= Vm*Vm; Vm4 = Vm2*Vm2;
753     Vm0     = PetscSqrtScalar(v0[2*lbus[i]]*v0[2*lbus[i]] + v0[2*lbus[i]+1]*v0[2*lbus[i]+1]);
754     PD      = QD = 0.0;
755     dPD_dVr = dPD_dVi = dQD_dVr = dQD_dVi = 0.0;
756     for (k=0; k < ld_nsegsp[i]; k++) {
757       PD      += ld_alphap[k]*PD0[i]*PetscPowScalar((Vm/Vm0),ld_betap[k]);
758       dPD_dVr += ld_alphap[k]*ld_betap[k]*PD0[i]*PetscPowScalar((1/Vm0),ld_betap[k])*Vr*PetscPowScalar(Vm,(ld_betap[k]-2));
759       dPD_dVi += ld_alphap[k]*ld_betap[k]*PD0[i]*PetscPowScalar((1/Vm0),ld_betap[k])*Vi*PetscPowScalar(Vm,(ld_betap[k]-2));
760     }
761     for (k=0; k < ld_nsegsq[i]; k++) {
762       QD      += ld_alphaq[k]*QD0[i]*PetscPowScalar((Vm/Vm0),ld_betaq[k]);
763       dQD_dVr += ld_alphaq[k]*ld_betaq[k]*QD0[i]*PetscPowScalar((1/Vm0),ld_betaq[k])*Vr*PetscPowScalar(Vm,(ld_betaq[k]-2));
764       dQD_dVi += ld_alphaq[k]*ld_betaq[k]*QD0[i]*PetscPowScalar((1/Vm0),ld_betaq[k])*Vi*PetscPowScalar(Vm,(ld_betaq[k]-2));
765     }
766 
767     /*    IDr = (PD*Vr + QD*Vi)/Vm2; */
768     /*    IDi = (-QD*Vr + PD*Vi)/Vm2; */
769 
770     dIDr_dVr = (dPD_dVr*Vr + dQD_dVr*Vi + PD)/Vm2 - ((PD*Vr + QD*Vi)*2*Vr)/Vm4;
771     dIDr_dVi = (dPD_dVi*Vr + dQD_dVi*Vi + QD)/Vm2 - ((PD*Vr + QD*Vi)*2*Vi)/Vm4;
772 
773     dIDi_dVr = (-dQD_dVr*Vr + dPD_dVr*Vi - QD)/Vm2 - ((-QD*Vr + PD*Vi)*2*Vr)/Vm4;
774     dIDi_dVi = (-dQD_dVi*Vr + dPD_dVi*Vi + PD)/Vm2 - ((-QD*Vr + PD*Vi)*2*Vi)/Vm4;
775 
776 
777     /*    fnet[2*lbus[i]]   += IDi; */
778     row[0] = net_start + 2*lbus[i];
779     col[0] = net_start + 2*lbus[i];  col[1] = net_start + 2*lbus[i]+1;
780     val[0] = dIDi_dVr;               val[1] = dIDi_dVi;
781     ierr   = MatSetValues(J,1,row,2,col,val,ADD_VALUES);CHKERRQ(ierr);
782     /*    fnet[2*lbus[i]+1] += IDr; */
783     row[0] = net_start + 2*lbus[i]+1;
784     col[0] = net_start + 2*lbus[i];  col[1] = net_start + 2*lbus[i]+1;
785     val[0] = dIDr_dVr;               val[1] = dIDr_dVi;
786     ierr   = MatSetValues(J,1,row,2,col,val,ADD_VALUES);CHKERRQ(ierr);
787   }
788   ierr = VecRestoreArray(user->V0,&v0);CHKERRQ(ierr);
789 
790   ierr = VecRestoreArray(Xgen,&xgen);CHKERRQ(ierr);
791   ierr = VecRestoreArray(Xnet,&xnet);CHKERRQ(ierr);
792 
793   ierr = DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);CHKERRQ(ierr);
794 
795   ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
796   ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
797   PetscFunctionReturn(0);
798 }
799 
800 /*
801    J = [I, 0
802         dg_dx, dg_dy]
803 */
804 PetscErrorCode AlgJacobian(SNES snes,Vec X,Mat A,Mat B,void *ctx)
805 {
806   PetscErrorCode ierr;
807   Userctx        *user=(Userctx*)ctx;
808 
809   PetscFunctionBegin;
810   ierr = ResidualJacobian(snes,X,A,B,ctx);CHKERRQ(ierr);
811   ierr = MatSetOption(A,MAT_KEEP_NONZERO_PATTERN,PETSC_TRUE);CHKERRQ(ierr);
812   ierr = MatZeroRowsIS(A,user->is_diff,1.0,NULL,NULL);CHKERRQ(ierr);
813   PetscFunctionReturn(0);
814 }
815 
816 /*
817    J = [a*I-df_dx, -df_dy
818         dg_dx, dg_dy]
819 */
820 
821 PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,Userctx *user)
822 {
823   PetscErrorCode ierr;
824   SNES           snes;
825   PetscScalar    atmp = (PetscScalar) a;
826   PetscInt       i,row;
827 
828   PetscFunctionBegin;
829   user->t = t;
830 
831   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
832   ierr = ResidualJacobian(snes,X,A,B,user);CHKERRQ(ierr);
833   for (i=0;i < ngen;i++) {
834     row = 9*i;
835     ierr = MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);CHKERRQ(ierr);
836     row  = 9*i+1;
837     ierr = MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);CHKERRQ(ierr);
838     row  = 9*i+2;
839     ierr = MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);CHKERRQ(ierr);
840     row  = 9*i+3;
841     ierr = MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);CHKERRQ(ierr);
842     row  = 9*i+6;
843     ierr = MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);CHKERRQ(ierr);
844     row  = 9*i+7;
845     ierr = MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);CHKERRQ(ierr);
846     row  = 9*i+8;
847     ierr = MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);CHKERRQ(ierr);
848   }
849   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
850   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
851   PetscFunctionReturn(0);
852 }
853 
854 /* Matrix JacobianP is constant so that it only needs to be evaluated once */
855 static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A, void *ctx0)
856 {
857   PetscErrorCode ierr;
858   PetscScalar    a;
859   PetscInt       row,col;
860   Userctx        *ctx=(Userctx*)ctx0;
861 
862   PetscFunctionBeginUser;
863 
864   if (ctx->jacp_flg) {
865     ierr = MatZeroEntries(A);CHKERRQ(ierr);
866 
867     for (col=0;col<3;col++) {
868       a    = 1.0/M[col];
869       row  = 9*col+3;
870       ierr = MatSetValues(A,1,&row,1,&col,&a,INSERT_VALUES);CHKERRQ(ierr);
871     }
872 
873     ctx->jacp_flg = PETSC_FALSE;
874 
875     ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
876     ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
877   }
878   PetscFunctionReturn(0);
879 }
880 
881 static PetscErrorCode CostIntegrand(TS ts,PetscReal t,Vec U,Vec R,Userctx *user)
882 {
883   PetscErrorCode    ierr;
884   const PetscScalar *u;
885   PetscInt          idx;
886   Vec               Xgen,Xnet;
887   PetscScalar       *r,*xgen;
888   PetscInt          i;
889 
890   PetscFunctionBegin;
891   ierr = DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);CHKERRQ(ierr);
892   ierr = DMCompositeScatter(user->dmpgrid,U,Xgen,Xnet);CHKERRQ(ierr);
893 
894   ierr = VecGetArray(Xgen,&xgen);CHKERRQ(ierr);
895 
896   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
897   ierr = VecGetArray(R,&r);CHKERRQ(ierr);
898   r[0] = 0.;
899   idx = 0;
900   for (i=0;i<ngen;i++) {
901     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);
902     idx  += 9;
903   }
904   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
905   ierr = VecRestoreArray(R,&r);CHKERRQ(ierr);
906   ierr = DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);CHKERRQ(ierr);
907   PetscFunctionReturn(0);
908 }
909 
910 static PetscErrorCode DRDUJacobianTranspose(TS ts,PetscReal t,Vec U,Mat DRDU,Mat B,Userctx *user)
911 {
912   PetscErrorCode ierr;
913   Vec            Xgen,Xnet,Dgen,Dnet;
914   PetscScalar    *xgen,*dgen;
915   PetscInt       i;
916   PetscInt       idx;
917   Vec            drdu_col;
918   PetscScalar    *xarr;
919 
920   PetscFunctionBegin;
921   ierr = VecDuplicate(U,&drdu_col);CHKERRQ(ierr);
922   ierr = MatDenseGetColumn(DRDU,0,&xarr);CHKERRQ(ierr);
923   ierr = VecPlaceArray(drdu_col,xarr);CHKERRQ(ierr);
924   ierr = DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);CHKERRQ(ierr);
925   ierr = DMCompositeGetLocalVectors(user->dmpgrid,&Dgen,&Dnet);CHKERRQ(ierr);
926   ierr = DMCompositeScatter(user->dmpgrid,U,Xgen,Xnet);CHKERRQ(ierr);
927   ierr = DMCompositeScatter(user->dmpgrid,drdu_col,Dgen,Dnet);CHKERRQ(ierr);
928 
929   ierr = VecGetArray(Xgen,&xgen);CHKERRQ(ierr);
930   ierr = VecGetArray(Dgen,&dgen);CHKERRQ(ierr);
931 
932   idx = 0;
933   for (i=0;i<ngen;i++) {
934     dgen[idx+3] = 0.;
935     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);
936     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);
937     idx += 9;
938   }
939 
940   ierr = VecRestoreArray(Dgen,&dgen);CHKERRQ(ierr);
941   ierr = VecRestoreArray(Xgen,&xgen);CHKERRQ(ierr);
942   ierr = DMCompositeGather(user->dmpgrid,INSERT_VALUES,drdu_col,Dgen,Dnet);CHKERRQ(ierr);
943   ierr = DMCompositeRestoreLocalVectors(user->dmpgrid,&Dgen,&Dnet);CHKERRQ(ierr);
944   ierr = DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);CHKERRQ(ierr);
945   ierr = VecResetArray(drdu_col);CHKERRQ(ierr);
946   ierr = MatDenseRestoreColumn(DRDU,&xarr);CHKERRQ(ierr);
947   ierr = VecDestroy(&drdu_col);CHKERRQ(ierr);
948   PetscFunctionReturn(0);
949 }
950 
951 static PetscErrorCode DRDPJacobianTranspose(TS ts,PetscReal t,Vec U,Mat drdp,Userctx *user)
952 {
953   PetscFunctionBegin;
954   PetscFunctionReturn(0);
955 }
956 
957 PetscErrorCode ComputeSensiP(Vec lambda,Vec mu,Vec *DICDP,Userctx *user)
958 {
959   PetscErrorCode ierr;
960   PetscScalar    *x,*y,sensip;
961   PetscInt       i;
962 
963   PetscFunctionBegin;
964   ierr = VecGetArray(lambda,&x);CHKERRQ(ierr);
965   ierr = VecGetArray(mu,&y);CHKERRQ(ierr);
966 
967   for (i=0;i<3;i++) {
968     ierr   = VecDot(lambda,DICDP[i],&sensip);CHKERRQ(ierr);
969     sensip = sensip+y[i];
970     /* ierr   = PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt %D th parameter: %g \n",i,(double)sensip);CHKERRQ(ierr); */
971      y[i] = sensip;
972   }
973   ierr = VecRestoreArray(mu,&y);CHKERRQ(ierr);
974   PetscFunctionReturn(0);
975 }
976 
977 int main(int argc,char **argv)
978 {
979   Userctx            user;
980   Vec                p;
981   PetscScalar        *x_ptr;
982   PetscErrorCode     ierr;
983   PetscMPIInt        size;
984   PetscInt           i;
985   PetscViewer        Xview,Ybusview;
986   PetscInt           *idx2;
987   Tao                tao;
988   KSP                ksp;
989   PC                 pc;
990   Vec                lowerb,upperb;
991 
992   ierr = PetscInitialize(&argc,&argv,"petscoptions",help);if (ierr) return ierr;
993   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
994   if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs");
995 
996   user.jacp_flg   = PETSC_TRUE;
997   user.neqs_gen   = 9*ngen; /* # eqs. for generator subsystem */
998   user.neqs_net   = 2*nbus; /* # eqs. for network subsystem   */
999   user.neqs_pgrid = user.neqs_gen + user.neqs_net;
1000 
1001   /* Create indices for differential and algebraic equations */
1002   ierr = PetscMalloc1(7*ngen,&idx2);CHKERRQ(ierr);
1003   for (i=0; i<ngen; i++) {
1004     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;
1005     idx2[7*i+4] = 9*i+6; idx2[7*i+5] = 9*i+7; idx2[7*i+6] = 9*i+8;
1006   }
1007   ierr = ISCreateGeneral(PETSC_COMM_WORLD,7*ngen,idx2,PETSC_COPY_VALUES,&user.is_diff);CHKERRQ(ierr);
1008   ierr = ISComplement(user.is_diff,0,user.neqs_pgrid,&user.is_alg);CHKERRQ(ierr);
1009   ierr = PetscFree(idx2);CHKERRQ(ierr);
1010 
1011   /* Set run time options */
1012   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Transient stability fault options","");CHKERRQ(ierr);
1013   {
1014     user.tfaulton  = 1.0;
1015     user.tfaultoff = 1.2;
1016     user.Rfault    = 0.0001;
1017     user.faultbus  = 8;
1018     ierr           = PetscOptionsReal("-tfaulton","","",user.tfaulton,&user.tfaulton,NULL);CHKERRQ(ierr);
1019     ierr           = PetscOptionsReal("-tfaultoff","","",user.tfaultoff,&user.tfaultoff,NULL);CHKERRQ(ierr);
1020     ierr           = PetscOptionsInt("-faultbus","","",user.faultbus,&user.faultbus,NULL);CHKERRQ(ierr);
1021     user.t0        = 0.0;
1022     user.tmax      = 1.3;
1023     ierr           = PetscOptionsReal("-t0","","",user.t0,&user.t0,NULL);CHKERRQ(ierr);
1024     ierr           = PetscOptionsReal("-tmax","","",user.tmax,&user.tmax,NULL);CHKERRQ(ierr);
1025     user.freq_u    = 61.0;
1026     user.freq_l    = 59.0;
1027     user.pow       = 2;
1028     ierr           = PetscOptionsReal("-frequ","","",user.freq_u,&user.freq_u,NULL);CHKERRQ(ierr);
1029     ierr           = PetscOptionsReal("-freql","","",user.freq_l,&user.freq_l,NULL);CHKERRQ(ierr);
1030     ierr           = PetscOptionsInt("-pow","","",user.pow,&user.pow,NULL);CHKERRQ(ierr);
1031 
1032   }
1033   ierr = PetscOptionsEnd();CHKERRQ(ierr);
1034 
1035   /* Create DMs for generator and network subsystems */
1036   ierr = DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,user.neqs_gen,1,1,NULL,&user.dmgen);CHKERRQ(ierr);
1037   ierr = DMSetOptionsPrefix(user.dmgen,"dmgen_");CHKERRQ(ierr);
1038   ierr = DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,user.neqs_net,1,1,NULL,&user.dmnet);CHKERRQ(ierr);
1039   ierr = DMSetOptionsPrefix(user.dmnet,"dmnet_");CHKERRQ(ierr);
1040   ierr = DMSetFromOptions(user.dmnet);CHKERRQ(ierr);
1041   ierr = DMSetUp(user.dmnet);CHKERRQ(ierr);
1042   /* Create a composite DM packer and add the two DMs */
1043   ierr = DMCompositeCreate(PETSC_COMM_WORLD,&user.dmpgrid);CHKERRQ(ierr);
1044   ierr = DMSetOptionsPrefix(user.dmpgrid,"pgrid_");CHKERRQ(ierr);
1045   ierr = DMSetFromOptions(user.dmgen);CHKERRQ(ierr);
1046   ierr = DMSetUp(user.dmgen);CHKERRQ(ierr);
1047   ierr = DMCompositeAddDM(user.dmpgrid,user.dmgen);CHKERRQ(ierr);
1048   ierr = DMCompositeAddDM(user.dmpgrid,user.dmnet);CHKERRQ(ierr);
1049 
1050   /* Read initial voltage vector and Ybus */
1051   ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,"X.bin",FILE_MODE_READ,&Xview);CHKERRQ(ierr);
1052   ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,"Ybus.bin",FILE_MODE_READ,&Ybusview);CHKERRQ(ierr);
1053 
1054   ierr = VecCreate(PETSC_COMM_WORLD,&user.V0);CHKERRQ(ierr);
1055   ierr = VecSetSizes(user.V0,PETSC_DECIDE,user.neqs_net);CHKERRQ(ierr);
1056   ierr = VecLoad(user.V0,Xview);CHKERRQ(ierr);
1057 
1058   ierr = MatCreate(PETSC_COMM_WORLD,&user.Ybus);CHKERRQ(ierr);
1059   ierr = MatSetSizes(user.Ybus,PETSC_DECIDE,PETSC_DECIDE,user.neqs_net,user.neqs_net);CHKERRQ(ierr);
1060   ierr = MatSetType(user.Ybus,MATBAIJ);CHKERRQ(ierr);
1061   /*  ierr = MatSetBlockSize(ctx->Ybus,2);CHKERRQ(ierr); */
1062   ierr = MatLoad(user.Ybus,Ybusview);CHKERRQ(ierr);
1063 
1064   ierr = PetscViewerDestroy(&Xview);CHKERRQ(ierr);
1065   ierr = PetscViewerDestroy(&Ybusview);CHKERRQ(ierr);
1066 
1067   /* Allocate space for Jacobians */
1068   ierr = MatCreate(PETSC_COMM_WORLD,&user.J);CHKERRQ(ierr);
1069   ierr = MatSetSizes(user.J,PETSC_DECIDE,PETSC_DECIDE,user.neqs_pgrid,user.neqs_pgrid);CHKERRQ(ierr);
1070   ierr = MatSetFromOptions(user.J);CHKERRQ(ierr);
1071   ierr = PreallocateJacobian(user.J,&user);CHKERRQ(ierr);
1072 
1073   ierr = MatCreate(PETSC_COMM_WORLD,&user.Jacp);CHKERRQ(ierr);
1074   ierr = MatSetSizes(user.Jacp,PETSC_DECIDE,PETSC_DECIDE,user.neqs_pgrid,3);CHKERRQ(ierr);
1075   ierr = MatSetFromOptions(user.Jacp);CHKERRQ(ierr);
1076   ierr = MatSetUp(user.Jacp);CHKERRQ(ierr);
1077   ierr = MatZeroEntries(user.Jacp);CHKERRQ(ierr); /* initialize to zeros */
1078 
1079   ierr = MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,3,1,NULL,&user.DRDP);CHKERRQ(ierr);
1080   ierr = MatSetUp(user.DRDP);CHKERRQ(ierr);
1081   ierr = MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,user.neqs_pgrid,1,NULL,&user.DRDU);CHKERRQ(ierr);
1082   ierr = MatSetUp(user.DRDU);CHKERRQ(ierr);
1083 
1084   /* Create TAO solver and set desired solution method */
1085   ierr = TaoCreate(PETSC_COMM_WORLD,&tao);CHKERRQ(ierr);
1086   ierr = TaoSetType(tao,TAOBLMVM);CHKERRQ(ierr);
1087   /*
1088      Optimization starts
1089   */
1090   /* Set initial solution guess */
1091   ierr = VecCreateSeq(PETSC_COMM_WORLD,3,&p);CHKERRQ(ierr);
1092   ierr = VecGetArray(p,&x_ptr);CHKERRQ(ierr);
1093   x_ptr[0] = PG[0]; x_ptr[1] = PG[1]; x_ptr[2] = PG[2];
1094   ierr = VecRestoreArray(p,&x_ptr);CHKERRQ(ierr);
1095 
1096   ierr = TaoSetInitialVector(tao,p);CHKERRQ(ierr);
1097   /* Set routine for function and gradient evaluation */
1098   ierr = TaoSetObjectiveAndGradientRoutine(tao,FormFunctionGradient,&user);CHKERRQ(ierr);
1099 
1100   /* Set bounds for the optimization */
1101   ierr = VecDuplicate(p,&lowerb);CHKERRQ(ierr);
1102   ierr = VecDuplicate(p,&upperb);CHKERRQ(ierr);
1103   ierr = VecGetArray(lowerb,&x_ptr);CHKERRQ(ierr);
1104   x_ptr[0] = 0.5; x_ptr[1] = 0.5; x_ptr[2] = 0.5;
1105   ierr = VecRestoreArray(lowerb,&x_ptr);CHKERRQ(ierr);
1106   ierr = VecGetArray(upperb,&x_ptr);CHKERRQ(ierr);
1107   x_ptr[0] = 2.0; x_ptr[1] = 2.0; x_ptr[2] = 2.0;
1108   ierr = VecRestoreArray(upperb,&x_ptr);CHKERRQ(ierr);
1109   ierr = TaoSetVariableBounds(tao,lowerb,upperb);CHKERRQ(ierr);
1110 
1111   /* Check for any TAO command line options */
1112   ierr = TaoSetFromOptions(tao);CHKERRQ(ierr);
1113   ierr = TaoGetKSP(tao,&ksp);CHKERRQ(ierr);
1114   if (ksp) {
1115     ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
1116     ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
1117   }
1118 
1119   /* SOLVE THE APPLICATION */
1120   ierr = TaoSolve(tao); CHKERRQ(ierr);
1121 
1122   ierr = VecView(p,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1123   /* Free TAO data structures */
1124   ierr = TaoDestroy(&tao);CHKERRQ(ierr);
1125 
1126   ierr = DMDestroy(&user.dmgen);CHKERRQ(ierr);
1127   ierr = DMDestroy(&user.dmnet);CHKERRQ(ierr);
1128   ierr = DMDestroy(&user.dmpgrid);CHKERRQ(ierr);
1129   ierr = ISDestroy(&user.is_diff);CHKERRQ(ierr);
1130   ierr = ISDestroy(&user.is_alg);CHKERRQ(ierr);
1131 
1132   ierr = MatDestroy(&user.J);CHKERRQ(ierr);
1133   ierr = MatDestroy(&user.Jacp);CHKERRQ(ierr);
1134   ierr = MatDestroy(&user.Ybus);CHKERRQ(ierr);
1135   /* ierr = MatDestroy(&user.Sol);CHKERRQ(ierr); */
1136   ierr = VecDestroy(&user.V0);CHKERRQ(ierr);
1137   ierr = VecDestroy(&p);CHKERRQ(ierr);
1138   ierr = VecDestroy(&lowerb);CHKERRQ(ierr);
1139   ierr = VecDestroy(&upperb);CHKERRQ(ierr);
1140   ierr = MatDestroy(&user.DRDU);CHKERRQ(ierr);
1141   ierr = MatDestroy(&user.DRDP);CHKERRQ(ierr);
1142   ierr = PetscFinalize();
1143   return ierr;
1144 }
1145 
1146 /* ------------------------------------------------------------------ */
1147 /*
1148    FormFunction - Evaluates the function and corresponding gradient.
1149 
1150    Input Parameters:
1151    tao - the Tao context
1152    X   - the input vector
1153    ptr - optional user-defined context, as set by TaoSetObjectiveAndGradientRoutine()
1154 
1155    Output Parameters:
1156    f   - the newly evaluated function
1157    G   - the newly evaluated gradient
1158 */
1159 PetscErrorCode FormFunctionGradient(Tao tao,Vec P,PetscReal *f,Vec G,void *ctx0)
1160 {
1161   TS             ts,quadts;
1162   SNES           snes_alg;
1163   PetscErrorCode ierr;
1164   Userctx        *ctx = (Userctx*)ctx0;
1165   Vec            X;
1166   PetscInt       i;
1167   /* sensitivity context */
1168   PetscScalar    *x_ptr;
1169   Vec            lambda[1],q;
1170   Vec            mu[1];
1171   PetscInt       steps1,steps2,steps3;
1172   Vec            DICDP[3];
1173   Vec            F_alg;
1174   PetscInt       row_loc,col_loc;
1175   PetscScalar    val;
1176   Vec            Xdot;
1177 
1178   PetscFunctionBegin;
1179   ierr  = VecGetArrayRead(P,(const PetscScalar**)&x_ptr);CHKERRQ(ierr);
1180   PG[0] = x_ptr[0];
1181   PG[1] = x_ptr[1];
1182   PG[2] = x_ptr[2];
1183   ierr  = VecRestoreArrayRead(P,(const PetscScalar**)&x_ptr);CHKERRQ(ierr);
1184 
1185   ctx->stepnum = 0;
1186 
1187   ierr = DMCreateGlobalVector(ctx->dmpgrid,&X);CHKERRQ(ierr);
1188 
1189   /* Create matrix to save solutions at each time step */
1190   /* ierr = MatCreateSeqDense(PETSC_COMM_SELF,ctx->neqs_pgrid+1,1002,NULL,&ctx->Sol);CHKERRQ(ierr); */
1191   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1192      Create timestepping solver context
1193      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1194   ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
1195   ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr);
1196   ierr = TSSetType(ts,TSCN);CHKERRQ(ierr);
1197   ierr = TSSetIFunction(ts,NULL,(TSIFunction) IFunction,ctx);CHKERRQ(ierr);
1198   ierr = TSSetIJacobian(ts,ctx->J,ctx->J,(TSIJacobian)IJacobian,ctx);CHKERRQ(ierr);
1199   ierr = TSSetApplicationContext(ts,ctx);CHKERRQ(ierr);
1200   /*   Set RHS JacobianP */
1201   ierr = TSSetRHSJacobianP(ts,ctx->Jacp,RHSJacobianP,ctx);CHKERRQ(ierr);
1202 
1203   ierr = TSCreateQuadratureTS(ts,PETSC_FALSE,&quadts);CHKERRQ(ierr);
1204   ierr = TSSetRHSFunction(quadts,NULL,(TSRHSFunction)CostIntegrand,ctx);CHKERRQ(ierr);
1205   ierr = TSSetRHSJacobian(quadts,ctx->DRDU,ctx->DRDU,(TSRHSJacobian)DRDUJacobianTranspose,ctx);CHKERRQ(ierr);
1206   ierr = TSSetRHSJacobianP(quadts,ctx->DRDP,(TSRHSJacobianP)DRDPJacobianTranspose,ctx);CHKERRQ(ierr);
1207 
1208   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1209      Set initial conditions
1210    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1211   ierr = SetInitialGuess(X,ctx);CHKERRQ(ierr);
1212 
1213   /* Approximate DICDP with finite difference, we want to zero out network variables */
1214   for (i=0;i<3;i++) {
1215     ierr = VecDuplicate(X,&DICDP[i]);CHKERRQ(ierr);
1216   }
1217   ierr = DICDPFiniteDifference(X,DICDP,ctx);CHKERRQ(ierr);
1218 
1219   ierr = VecDuplicate(X,&F_alg);CHKERRQ(ierr);
1220   ierr = SNESCreate(PETSC_COMM_WORLD,&snes_alg);CHKERRQ(ierr);
1221   ierr = SNESSetFunction(snes_alg,F_alg,AlgFunction,ctx);CHKERRQ(ierr);
1222   ierr = MatZeroEntries(ctx->J);CHKERRQ(ierr);
1223   ierr = SNESSetJacobian(snes_alg,ctx->J,ctx->J,AlgJacobian,ctx);CHKERRQ(ierr);
1224   ierr = SNESSetOptionsPrefix(snes_alg,"alg_");CHKERRQ(ierr);
1225   ierr = SNESSetFromOptions(snes_alg);CHKERRQ(ierr);
1226   ctx->alg_flg = PETSC_TRUE;
1227   /* Solve the algebraic equations */
1228   ierr = SNESSolve(snes_alg,NULL,X);CHKERRQ(ierr);
1229 
1230   /* Just to set up the Jacobian structure */
1231   ierr = VecDuplicate(X,&Xdot);CHKERRQ(ierr);
1232   ierr = IJacobian(ts,0.0,X,Xdot,0.0,ctx->J,ctx->J,ctx);CHKERRQ(ierr);
1233   ierr = VecDestroy(&Xdot);CHKERRQ(ierr);
1234 
1235   ctx->stepnum++;
1236 
1237   /*
1238     Save trajectory of solution so that TSAdjointSolve() may be used
1239   */
1240   ierr = TSSetSaveTrajectory(ts);CHKERRQ(ierr);
1241 
1242   ierr = TSSetTimeStep(ts,0.01);CHKERRQ(ierr);
1243   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr);
1244   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
1245   /* ierr = TSSetPostStep(ts,SaveSolution);CHKERRQ(ierr); */
1246 
1247 
1248   /* Prefault period */
1249   ctx->alg_flg = PETSC_FALSE;
1250   ierr = TSSetTime(ts,0.0);CHKERRQ(ierr);
1251   ierr = TSSetMaxTime(ts,ctx->tfaulton);CHKERRQ(ierr);
1252   ierr = TSSolve(ts,X);CHKERRQ(ierr);
1253   ierr = TSGetStepNumber(ts,&steps1);CHKERRQ(ierr);
1254 
1255   /* Create the nonlinear solver for solving the algebraic system */
1256   /* Note that although the algebraic system needs to be solved only for
1257      Idq and V, we reuse the entire system including xgen. The xgen
1258      variables are held constant by setting their residuals to 0 and
1259      putting a 1 on the Jacobian diagonal for xgen rows
1260   */
1261   ierr = MatZeroEntries(ctx->J);CHKERRQ(ierr);
1262 
1263   /* Apply disturbance - resistive fault at ctx->faultbus */
1264   /* This is done by adding shunt conductance to the diagonal location
1265      in the Ybus matrix */
1266   row_loc = 2*ctx->faultbus; col_loc = 2*ctx->faultbus+1; /* Location for G */
1267   val     = 1/ctx->Rfault;
1268   ierr    = MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);CHKERRQ(ierr);
1269   row_loc = 2*ctx->faultbus+1; col_loc = 2*ctx->faultbus; /* Location for G */
1270   val     = 1/ctx->Rfault;
1271   ierr    = MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);CHKERRQ(ierr);
1272 
1273   ierr = MatAssemblyBegin(ctx->Ybus,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1274   ierr = MatAssemblyEnd(ctx->Ybus,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1275 
1276   ctx->alg_flg = PETSC_TRUE;
1277   /* Solve the algebraic equations */
1278   ierr = SNESSolve(snes_alg,NULL,X);CHKERRQ(ierr);
1279 
1280   ctx->stepnum++;
1281 
1282   /* Disturbance period */
1283   ctx->alg_flg = PETSC_FALSE;
1284   ierr = TSSetTime(ts,ctx->tfaulton);CHKERRQ(ierr);
1285   ierr = TSSetMaxTime(ts,ctx->tfaultoff);CHKERRQ(ierr);
1286   ierr = TSSolve(ts,X);CHKERRQ(ierr);
1287   ierr = TSGetStepNumber(ts,&steps2);CHKERRQ(ierr);
1288   steps2 -= steps1;
1289 
1290   /* Remove the fault */
1291   row_loc = 2*ctx->faultbus; col_loc = 2*ctx->faultbus+1;
1292   val     = -1/ctx->Rfault;
1293   ierr    = MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);CHKERRQ(ierr);
1294   row_loc = 2*ctx->faultbus+1; col_loc = 2*ctx->faultbus;
1295   val     = -1/ctx->Rfault;
1296   ierr    = MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);CHKERRQ(ierr);
1297 
1298   ierr = MatAssemblyBegin(ctx->Ybus,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1299   ierr = MatAssemblyEnd(ctx->Ybus,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1300 
1301   ierr = MatZeroEntries(ctx->J);CHKERRQ(ierr);
1302 
1303   ctx->alg_flg = PETSC_TRUE;
1304 
1305   /* Solve the algebraic equations */
1306   ierr = SNESSolve(snes_alg,NULL,X);CHKERRQ(ierr);
1307 
1308   ctx->stepnum++;
1309 
1310   /* Post-disturbance period */
1311   ctx->alg_flg = PETSC_TRUE;
1312   ierr = TSSetTime(ts,ctx->tfaultoff);CHKERRQ(ierr);
1313   ierr = TSSetMaxTime(ts,ctx->tmax);CHKERRQ(ierr);
1314   ierr = TSSolve(ts,X);CHKERRQ(ierr);
1315   ierr = TSGetStepNumber(ts,&steps3);CHKERRQ(ierr);
1316   steps3 -= steps2;
1317   steps3 -= steps1;
1318 
1319   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1320      Adjoint model starts here
1321      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1322   ierr = TSSetPostStep(ts,NULL);CHKERRQ(ierr);
1323   ierr = MatCreateVecs(ctx->J,&lambda[0],NULL);CHKERRQ(ierr);
1324   /*   Set initial conditions for the adjoint integration */
1325   ierr = VecZeroEntries(lambda[0]);CHKERRQ(ierr);
1326 
1327   ierr = MatCreateVecs(ctx->Jacp,&mu[0],NULL);CHKERRQ(ierr);
1328   ierr = VecZeroEntries(mu[0]);CHKERRQ(ierr);
1329   ierr = TSSetCostGradients(ts,1,lambda,mu);CHKERRQ(ierr);
1330 
1331   ierr = TSAdjointSetSteps(ts,steps3);CHKERRQ(ierr);
1332   ierr = TSAdjointSolve(ts);CHKERRQ(ierr);
1333 
1334   ierr = MatZeroEntries(ctx->J);CHKERRQ(ierr);
1335   /* Applying disturbance - resistive fault at ctx->faultbus */
1336   /* This is done by deducting shunt conductance to the diagonal location
1337      in the Ybus matrix */
1338   row_loc = 2*ctx->faultbus; col_loc = 2*ctx->faultbus+1; /* Location for G */
1339   val     = 1./ctx->Rfault;
1340   ierr    = MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);CHKERRQ(ierr);
1341   row_loc = 2*ctx->faultbus+1; col_loc = 2*ctx->faultbus; /* Location for G */
1342   val     = 1./ctx->Rfault;
1343   ierr    = MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);CHKERRQ(ierr);
1344 
1345   ierr = MatAssemblyBegin(ctx->Ybus,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1346   ierr = MatAssemblyEnd(ctx->Ybus,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1347 
1348 
1349   /*   Set number of steps for the adjoint integration */
1350   ierr = TSAdjointSetSteps(ts,steps2);CHKERRQ(ierr);
1351   ierr = TSAdjointSolve(ts);CHKERRQ(ierr);
1352 
1353   ierr = MatZeroEntries(ctx->J);CHKERRQ(ierr);
1354   /* remove the fault */
1355   row_loc = 2*ctx->faultbus; col_loc = 2*ctx->faultbus+1; /* Location for G */
1356   val     = -1./ctx->Rfault;
1357   ierr    = MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);CHKERRQ(ierr);
1358   row_loc = 2*ctx->faultbus+1; col_loc = 2*ctx->faultbus; /* Location for G */
1359   val     = -1./ctx->Rfault;
1360   ierr    = MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);CHKERRQ(ierr);
1361 
1362   ierr = MatAssemblyBegin(ctx->Ybus,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1363   ierr = MatAssemblyEnd(ctx->Ybus,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1364 
1365   /*   Set number of steps for the adjoint integration */
1366   ierr = TSAdjointSetSteps(ts,steps1);CHKERRQ(ierr);
1367   ierr = TSAdjointSolve(ts);CHKERRQ(ierr);
1368 
1369   ierr = ComputeSensiP(lambda[0],mu[0],DICDP,ctx);CHKERRQ(ierr);
1370   ierr = VecCopy(mu[0],G);CHKERRQ(ierr);
1371 
1372   ierr = TSGetQuadratureTS(ts,NULL,&quadts);CHKERRQ(ierr);
1373   ierr = TSGetSolution(quadts,&q);CHKERRQ(ierr);
1374   ierr = VecGetArray(q,&x_ptr);CHKERRQ(ierr);
1375   *f   = x_ptr[0];
1376   x_ptr[0] = 0;
1377   ierr = VecRestoreArray(q,&x_ptr);CHKERRQ(ierr);
1378 
1379   ierr = VecDestroy(&lambda[0]);CHKERRQ(ierr);
1380   ierr = VecDestroy(&mu[0]);CHKERRQ(ierr);
1381 
1382   ierr = SNESDestroy(&snes_alg);CHKERRQ(ierr);
1383   ierr = VecDestroy(&F_alg);CHKERRQ(ierr);
1384   ierr = VecDestroy(&X);CHKERRQ(ierr);
1385   ierr = TSDestroy(&ts);CHKERRQ(ierr);
1386   for (i=0;i<3;i++) {
1387     ierr = VecDestroy(&DICDP[i]);CHKERRQ(ierr);
1388   }
1389   PetscFunctionReturn(0);
1390 }
1391 
1392 /*TEST
1393 
1394    build:
1395       requires: double !complex !define(PETSC_USE_64BIT_INDICES)
1396 
1397    test:
1398       args: -viewer_binary_skip_info -tao_monitor -tao_gttol .2
1399       localrunfiles: petscoptions X.bin Ybus.bin
1400 
1401 TEST*/
1402