xref: /petsc/src/ts/tutorials/power_grid/stability_9bus/ex9busoptfd.c (revision 8fb5bd83c3955fefcf33a54e3bb66920a9fa884b)
1 static char help[] = "Using finite difference for the problem in ex9busopt.c \n\n";
2 
3 /*
4   Use finite difference approximations to solve the same optimization problem as in ex9busopt.c.
5  */
6 
7 #include <petsctao.h>
8 #include <petscts.h>
9 #include <petscdm.h>
10 #include <petscdmda.h>
11 #include <petscdmcomposite.h>
12 
13 PetscErrorCode FormFunction(Tao,Vec,PetscReal*,void*);
14 
15 #define freq 60
16 #define w_s (2*PETSC_PI*freq)
17 
18 /* Sizes and indices */
19 const PetscInt nbus    = 9; /* Number of network buses */
20 const PetscInt ngen    = 3; /* Number of generators */
21 const PetscInt nload   = 3; /* Number of loads */
22 const PetscInt gbus[3] = {0,1,2}; /* Buses at which generators are incident */
23 const PetscInt lbus[3] = {4,5,7}; /* Buses at which loads are incident */
24 
25 /* Generator real and reactive powers (found via loadflow) */
26 PetscScalar PG[3] = { 0.69,1.59,0.69};
27 /* PetscScalar PG[3] = {0.716786142395021,1.630000000000000,0.850000000000000};*/
28 const PetscScalar QG[3] = {0.270702180178785,0.066120127797275,-0.108402221791588};
29 /* Generator constants */
30 const PetscScalar H[3]    = {23.64,6.4,3.01};   /* Inertia constant */
31 const PetscScalar Rs[3]   = {0.0,0.0,0.0}; /* Stator Resistance */
32 const PetscScalar Xd[3]   = {0.146,0.8958,1.3125};  /* d-axis reactance */
33 const PetscScalar Xdp[3]  = {0.0608,0.1198,0.1813}; /* d-axis transient reactance */
34 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 */
35 const PetscScalar Xqp[3]  = {0.0969,0.1969,0.25}; /* q-axis transient reactance */
36 const PetscScalar Td0p[3] = {8.96,6.0,5.89}; /* d-axis open circuit time constant */
37 const PetscScalar Tq0p[3] = {0.31,0.535,0.6}; /* q-axis open circuit time constant */
38 PetscScalar M[3]; /* M = 2*H/w_s */
39 PetscScalar D[3]; /* D = 0.1*M */
40 
41 PetscScalar TM[3]; /* Mechanical Torque */
42 /* Exciter system constants */
43 const PetscScalar KA[3] = {20.0,20.0,20.0};  /* Voltage regulartor gain constant */
44 const PetscScalar TA[3] = {0.2,0.2,0.2};     /* Voltage regulator time constant */
45 const PetscScalar KE[3] = {1.0,1.0,1.0};     /* Exciter gain constant */
46 const PetscScalar TE[3] = {0.314,0.314,0.314}; /* Exciter time constant */
47 const PetscScalar KF[3] = {0.063,0.063,0.063};  /* Feedback stabilizer gain constant */
48 const PetscScalar TF[3] = {0.35,0.35,0.35};    /* Feedback stabilizer time constant */
49 const PetscScalar k1[3] = {0.0039,0.0039,0.0039};
50 const PetscScalar k2[3] = {1.555,1.555,1.555};  /* k1 and k2 for calculating the saturation function SE = k1*exp(k2*Efd) */
51 
52 PetscScalar Vref[3];
53 /* Load constants
54   We use a composite load model that describes the load and reactive powers at each time instant as follows
55   P(t) = \sum\limits_{i=0}^ld_nsegsp \ld_alphap_i*P_D0(\frac{V_m(t)}{V_m0})^\ld_betap_i
56   Q(t) = \sum\limits_{i=0}^ld_nsegsq \ld_alphaq_i*Q_D0(\frac{V_m(t)}{V_m0})^\ld_betaq_i
57   where
58     ld_nsegsp,ld_nsegsq - Number of individual load models for real and reactive power loads
59     ld_alphap,ld_alphap - Percentage contribution (weights) or loads
60     P_D0                - Real power load
61     Q_D0                - Reactive power load
62     V_m(t)              - Voltage magnitude at time t
63     V_m0                - Voltage magnitude at t = 0
64     ld_betap, ld_betaq  - exponents describing the load model for real and reactive part
65 
66     Note: All loads have the same characteristic currently.
67 */
68 const PetscScalar PD0[3] = {1.25,0.9,1.0};
69 const PetscScalar QD0[3] = {0.5,0.3,0.35};
70 const PetscInt    ld_nsegsp[3] = {3,3,3};
71 const PetscScalar ld_alphap[3] = {1.0,0.0,0.0};
72 const PetscScalar ld_betap[3]  = {2.0,1.0,0.0};
73 const PetscInt    ld_nsegsq[3] = {3,3,3};
74 const PetscScalar ld_alphaq[3] = {1.0,0.0,0.0};
75 const PetscScalar ld_betaq[3]  = {2.0,1.0,0.0};
76 
77 typedef struct {
78   DM          dmgen, dmnet; /* DMs to manage generator and network subsystem */
79   DM          dmpgrid; /* Composite DM to manage the entire power grid */
80   Mat         Ybus; /* Network admittance matrix */
81   Vec         V0;  /* Initial voltage vector (Power flow solution) */
82   PetscReal   tfaulton,tfaultoff; /* Fault on and off times */
83   PetscInt    faultbus; /* Fault bus */
84   PetscScalar Rfault;
85   PetscReal   t0,tmax;
86   PetscInt    neqs_gen,neqs_net,neqs_pgrid;
87   Mat         Sol; /* Matrix to save solution at each time step */
88   PetscInt    stepnum;
89   PetscBool   alg_flg;
90   PetscReal   t;
91   IS          is_diff; /* indices for differential equations */
92   IS          is_alg; /* indices for algebraic equations */
93   PetscReal   freq_u,freq_l; /* upper and lower frequency limit */
94   PetscInt    pow; /* power coefficient used in the cost function */
95   Vec         vec_q;
96 } Userctx;
97 
98 /* Converts from machine frame (dq) to network (phase a real,imag) reference frame */
99 PetscErrorCode dq2ri(PetscScalar Fd,PetscScalar Fq,PetscScalar delta,PetscScalar *Fr, PetscScalar *Fi)
100 {
101   PetscFunctionBegin;
102   *Fr =  Fd*PetscSinScalar(delta) + Fq*PetscCosScalar(delta);
103   *Fi = -Fd*PetscCosScalar(delta) + Fq*PetscSinScalar(delta);
104   PetscFunctionReturn(0);
105 }
106 
107 /* Converts from network frame ([phase a real,imag) to machine (dq) reference frame */
108 PetscErrorCode ri2dq(PetscScalar Fr,PetscScalar Fi,PetscScalar delta,PetscScalar *Fd, PetscScalar *Fq)
109 {
110   PetscFunctionBegin;
111   *Fd =  Fr*PetscSinScalar(delta) - Fi*PetscCosScalar(delta);
112   *Fq =  Fr*PetscCosScalar(delta) + Fi*PetscSinScalar(delta);
113   PetscFunctionReturn(0);
114 }
115 
116 PetscErrorCode SetInitialGuess(Vec X,Userctx *user)
117 {
118   Vec            Xgen,Xnet;
119   PetscScalar    *xgen,*xnet;
120   PetscInt       i,idx=0;
121   PetscScalar    Vr,Vi,IGr,IGi,Vm,Vm2;
122   PetscScalar    Eqp,Edp,delta;
123   PetscScalar    Efd,RF,VR; /* Exciter variables */
124   PetscScalar    Id,Iq;  /* Generator dq axis currents */
125   PetscScalar    theta,Vd,Vq,SE;
126 
127   PetscFunctionBegin;
128   M[0] = 2*H[0]/w_s; M[1] = 2*H[1]/w_s; M[2] = 2*H[2]/w_s;
129   D[0] = 0.1*M[0]; D[1] = 0.1*M[1]; D[2] = 0.1*M[2];
130 
131   PetscCall(DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet));
132 
133   /* Network subsystem initialization */
134   PetscCall(VecCopy(user->V0,Xnet));
135 
136   /* Generator subsystem initialization */
137   PetscCall(VecGetArray(Xgen,&xgen));
138   PetscCall(VecGetArray(Xnet,&xnet));
139 
140   for (i=0; i < ngen; i++) {
141     Vr  = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */
142     Vi  = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */
143     Vm  = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm;
144     IGr = (Vr*PG[i] + Vi*QG[i])/Vm2;
145     IGi = (Vi*PG[i] - Vr*QG[i])/Vm2;
146 
147     delta = PetscAtan2Real(Vi+Xq[i]*IGr,Vr-Xq[i]*IGi); /* Machine angle */
148 
149     theta = PETSC_PI/2.0 - delta;
150 
151     Id = IGr*PetscCosScalar(theta) - IGi*PetscSinScalar(theta); /* d-axis stator current */
152     Iq = IGr*PetscSinScalar(theta) + IGi*PetscCosScalar(theta); /* q-axis stator current */
153 
154     Vd = Vr*PetscCosScalar(theta) - Vi*PetscSinScalar(theta);
155     Vq = Vr*PetscSinScalar(theta) + Vi*PetscCosScalar(theta);
156 
157     Edp = Vd + Rs[i]*Id - Xqp[i]*Iq; /* d-axis transient EMF */
158     Eqp = Vq + Rs[i]*Iq + Xdp[i]*Id; /* q-axis transient EMF */
159 
160     TM[i] = PG[i];
161 
162     /* The generator variables are ordered as [Eqp,Edp,delta,w,Id,Iq] */
163     xgen[idx]   = Eqp;
164     xgen[idx+1] = Edp;
165     xgen[idx+2] = delta;
166     xgen[idx+3] = w_s;
167 
168     idx = idx + 4;
169 
170     xgen[idx]   = Id;
171     xgen[idx+1] = Iq;
172 
173     idx = idx + 2;
174 
175     /* Exciter */
176     Efd = Eqp + (Xd[i] - Xdp[i])*Id;
177     SE  = k1[i]*PetscExpScalar(k2[i]*Efd);
178     VR  =  KE[i]*Efd + SE;
179     RF  =  KF[i]*Efd/TF[i];
180 
181     xgen[idx]   = Efd;
182     xgen[idx+1] = RF;
183     xgen[idx+2] = VR;
184 
185     Vref[i] = Vm + (VR/KA[i]);
186 
187     idx = idx + 3;
188   }
189 
190   PetscCall(VecRestoreArray(Xgen,&xgen));
191   PetscCall(VecRestoreArray(Xnet,&xnet));
192 
193   /* PetscCall(VecView(Xgen,0)); */
194   PetscCall(DMCompositeGather(user->dmpgrid,INSERT_VALUES,X,Xgen,Xnet));
195   PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet));
196   PetscFunctionReturn(0);
197 }
198 
199 /* Computes F = [-f(x,y);g(x,y)] */
200 PetscErrorCode ResidualFunction(SNES snes,Vec X, Vec F, Userctx *user)
201 {
202   Vec            Xgen,Xnet,Fgen,Fnet;
203   PetscScalar    *xgen,*xnet,*fgen,*fnet;
204   PetscInt       i,idx=0;
205   PetscScalar    Vr,Vi,Vm,Vm2;
206   PetscScalar    Eqp,Edp,delta,w; /* Generator variables */
207   PetscScalar    Efd,RF,VR; /* Exciter variables */
208   PetscScalar    Id,Iq;  /* Generator dq axis currents */
209   PetscScalar    Vd,Vq,SE;
210   PetscScalar    IGr,IGi,IDr,IDi;
211   PetscScalar    Zdq_inv[4],det;
212   PetscScalar    PD,QD,Vm0,*v0;
213   PetscInt       k;
214 
215   PetscFunctionBegin;
216   PetscCall(VecZeroEntries(F));
217   PetscCall(DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet));
218   PetscCall(DMCompositeGetLocalVectors(user->dmpgrid,&Fgen,&Fnet));
219   PetscCall(DMCompositeScatter(user->dmpgrid,X,Xgen,Xnet));
220   PetscCall(DMCompositeScatter(user->dmpgrid,F,Fgen,Fnet));
221 
222   /* Network current balance residual IG + Y*V + IL = 0. Only YV is added here.
223      The generator current injection, IG, and load current injection, ID are added later
224   */
225   /* Note that the values in Ybus are stored assuming the imaginary current balance
226      equation is ordered first followed by real current balance equation for each bus.
227      Thus imaginary current contribution goes in location 2*i, and
228      real current contribution in 2*i+1
229   */
230   PetscCall(MatMult(user->Ybus,Xnet,Fnet));
231 
232   PetscCall(VecGetArray(Xgen,&xgen));
233   PetscCall(VecGetArray(Xnet,&xnet));
234   PetscCall(VecGetArray(Fgen,&fgen));
235   PetscCall(VecGetArray(Fnet,&fnet));
236 
237   /* Generator subsystem */
238   for (i=0; i < ngen; i++) {
239     Eqp   = xgen[idx];
240     Edp   = xgen[idx+1];
241     delta = xgen[idx+2];
242     w     = xgen[idx+3];
243     Id    = xgen[idx+4];
244     Iq    = xgen[idx+5];
245     Efd   = xgen[idx+6];
246     RF    = xgen[idx+7];
247     VR    = xgen[idx+8];
248 
249     /* Generator differential equations */
250     fgen[idx]   = (Eqp + (Xd[i] - Xdp[i])*Id - Efd)/Td0p[i];
251     fgen[idx+1] = (Edp - (Xq[i] - Xqp[i])*Iq)/Tq0p[i];
252     fgen[idx+2] = -w + w_s;
253     fgen[idx+3] = (-TM[i] + Edp*Id + Eqp*Iq + (Xqp[i] - Xdp[i])*Id*Iq + D[i]*(w - w_s))/M[i];
254 
255     Vr = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */
256     Vi = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */
257 
258     PetscCall(ri2dq(Vr,Vi,delta,&Vd,&Vq));
259     /* Algebraic equations for stator currents */
260     det = Rs[i]*Rs[i] + Xdp[i]*Xqp[i];
261 
262     Zdq_inv[0] = Rs[i]/det;
263     Zdq_inv[1] = Xqp[i]/det;
264     Zdq_inv[2] = -Xdp[i]/det;
265     Zdq_inv[3] = Rs[i]/det;
266 
267     fgen[idx+4] = Zdq_inv[0]*(-Edp + Vd) + Zdq_inv[1]*(-Eqp + Vq) + Id;
268     fgen[idx+5] = Zdq_inv[2]*(-Edp + Vd) + Zdq_inv[3]*(-Eqp + Vq) + Iq;
269 
270     /* Add generator current injection to network */
271     PetscCall(dq2ri(Id,Iq,delta,&IGr,&IGi));
272 
273     fnet[2*gbus[i]]   -= IGi;
274     fnet[2*gbus[i]+1] -= IGr;
275 
276     Vm = PetscSqrtScalar(Vd*Vd + Vq*Vq);
277 
278     SE = k1[i]*PetscExpScalar(k2[i]*Efd);
279 
280     /* Exciter differential equations */
281     fgen[idx+6] = (KE[i]*Efd + SE - VR)/TE[i];
282     fgen[idx+7] = (RF - KF[i]*Efd/TF[i])/TF[i];
283     fgen[idx+8] = (VR - KA[i]*RF + KA[i]*KF[i]*Efd/TF[i] - KA[i]*(Vref[i] - Vm))/TA[i];
284 
285     idx = idx + 9;
286   }
287 
288   PetscCall(VecGetArray(user->V0,&v0));
289   for (i=0; i < nload; i++) {
290     Vr  = xnet[2*lbus[i]]; /* Real part of load bus voltage */
291     Vi  = xnet[2*lbus[i]+1]; /* Imaginary part of the load bus voltage */
292     Vm  = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm;
293     Vm0 = PetscSqrtScalar(v0[2*lbus[i]]*v0[2*lbus[i]] + v0[2*lbus[i]+1]*v0[2*lbus[i]+1]);
294     PD  = QD = 0.0;
295     for (k=0; k < ld_nsegsp[i]; k++) PD += ld_alphap[k]*PD0[i]*PetscPowScalar((Vm/Vm0),ld_betap[k]);
296     for (k=0; k < ld_nsegsq[i]; k++) QD += ld_alphaq[k]*QD0[i]*PetscPowScalar((Vm/Vm0),ld_betaq[k]);
297 
298     /* Load currents */
299     IDr = (PD*Vr + QD*Vi)/Vm2;
300     IDi = (-QD*Vr + PD*Vi)/Vm2;
301 
302     fnet[2*lbus[i]]   += IDi;
303     fnet[2*lbus[i]+1] += IDr;
304   }
305   PetscCall(VecRestoreArray(user->V0,&v0));
306 
307   PetscCall(VecRestoreArray(Xgen,&xgen));
308   PetscCall(VecRestoreArray(Xnet,&xnet));
309   PetscCall(VecRestoreArray(Fgen,&fgen));
310   PetscCall(VecRestoreArray(Fnet,&fnet));
311 
312   PetscCall(DMCompositeGather(user->dmpgrid,INSERT_VALUES,F,Fgen,Fnet));
313   PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet));
314   PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid,&Fgen,&Fnet));
315   PetscFunctionReturn(0);
316 }
317 
318 /* \dot{x} - f(x,y)
319      g(x,y) = 0
320  */
321 PetscErrorCode IFunction(TS ts,PetscReal t, Vec X, Vec Xdot, Vec F, Userctx *user)
322 {
323   SNES              snes;
324   PetscScalar       *f;
325   const PetscScalar *xdot;
326   PetscInt          i;
327 
328   PetscFunctionBegin;
329   user->t = t;
330 
331   PetscCall(TSGetSNES(ts,&snes));
332   PetscCall(ResidualFunction(snes,X,F,user));
333   PetscCall(VecGetArray(F,&f));
334   PetscCall(VecGetArrayRead(Xdot,&xdot));
335   for (i=0;i < ngen;i++) {
336     f[9*i]   += xdot[9*i];
337     f[9*i+1] += xdot[9*i+1];
338     f[9*i+2] += xdot[9*i+2];
339     f[9*i+3] += xdot[9*i+3];
340     f[9*i+6] += xdot[9*i+6];
341     f[9*i+7] += xdot[9*i+7];
342     f[9*i+8] += xdot[9*i+8];
343   }
344   PetscCall(VecRestoreArray(F,&f));
345   PetscCall(VecRestoreArrayRead(Xdot,&xdot));
346   PetscFunctionReturn(0);
347 }
348 
349 /* This function is used for solving the algebraic system only during fault on and
350    off times. It computes the entire F and then zeros out the part corresponding to
351    differential equations
352  F = [0;g(y)];
353 */
354 PetscErrorCode AlgFunction(SNES snes, Vec X, Vec F, void *ctx)
355 {
356   Userctx        *user=(Userctx*)ctx;
357   PetscScalar    *f;
358   PetscInt       i;
359 
360   PetscFunctionBegin;
361   PetscCall(ResidualFunction(snes,X,F,user));
362   PetscCall(VecGetArray(F,&f));
363   for (i=0; i < ngen; i++) {
364     f[9*i]   = 0;
365     f[9*i+1] = 0;
366     f[9*i+2] = 0;
367     f[9*i+3] = 0;
368     f[9*i+6] = 0;
369     f[9*i+7] = 0;
370     f[9*i+8] = 0;
371   }
372   PetscCall(VecRestoreArray(F,&f));
373   PetscFunctionReturn(0);
374 }
375 
376 PetscErrorCode PreallocateJacobian(Mat J, Userctx *user)
377 {
378   PetscInt       *d_nnz;
379   PetscInt       i,idx=0,start=0;
380   PetscInt       ncols;
381 
382   PetscFunctionBegin;
383   PetscCall(PetscMalloc1(user->neqs_pgrid,&d_nnz));
384   for (i=0; i<user->neqs_pgrid; i++) d_nnz[i] = 0;
385   /* Generator subsystem */
386   for (i=0; i < ngen; i++) {
387 
388     d_nnz[idx]   += 3;
389     d_nnz[idx+1] += 2;
390     d_nnz[idx+2] += 2;
391     d_nnz[idx+3] += 5;
392     d_nnz[idx+4] += 6;
393     d_nnz[idx+5] += 6;
394 
395     d_nnz[user->neqs_gen+2*gbus[i]]   += 3;
396     d_nnz[user->neqs_gen+2*gbus[i]+1] += 3;
397 
398     d_nnz[idx+6] += 2;
399     d_nnz[idx+7] += 2;
400     d_nnz[idx+8] += 5;
401 
402     idx = idx + 9;
403   }
404 
405   start = user->neqs_gen;
406 
407   for (i=0; i < nbus; i++) {
408     PetscCall(MatGetRow(user->Ybus,2*i,&ncols,NULL,NULL));
409     d_nnz[start+2*i]   += ncols;
410     d_nnz[start+2*i+1] += ncols;
411     PetscCall(MatRestoreRow(user->Ybus,2*i,&ncols,NULL,NULL));
412   }
413 
414   PetscCall(MatSeqAIJSetPreallocation(J,0,d_nnz));
415 
416   PetscCall(PetscFree(d_nnz));
417   PetscFunctionReturn(0);
418 }
419 
420 /*
421    J = [-df_dx, -df_dy
422         dg_dx, dg_dy]
423 */
424 PetscErrorCode ResidualJacobian(SNES snes,Vec X,Mat J,Mat B,void *ctx)
425 {
426   Userctx           *user=(Userctx*)ctx;
427   Vec               Xgen,Xnet;
428   PetscScalar       *xgen,*xnet;
429   PetscInt          i,idx=0;
430   PetscScalar       Vr,Vi,Vm,Vm2;
431   PetscScalar       Eqp,Edp,delta; /* Generator variables */
432   PetscScalar       Efd; /* Exciter variables */
433   PetscScalar       Id,Iq;  /* Generator dq axis currents */
434   PetscScalar       Vd,Vq;
435   PetscScalar       val[10];
436   PetscInt          row[2],col[10];
437   PetscInt          net_start=user->neqs_gen;
438   PetscScalar       Zdq_inv[4],det;
439   PetscScalar       dVd_dVr,dVd_dVi,dVq_dVr,dVq_dVi,dVd_ddelta,dVq_ddelta;
440   PetscScalar       dIGr_ddelta,dIGi_ddelta,dIGr_dId,dIGr_dIq,dIGi_dId,dIGi_dIq;
441   PetscScalar       dSE_dEfd;
442   PetscScalar       dVm_dVd,dVm_dVq,dVm_dVr,dVm_dVi;
443   PetscInt          ncols;
444   const PetscInt    *cols;
445   const PetscScalar *yvals;
446   PetscInt          k;
447   PetscScalar       PD,QD,Vm0,*v0,Vm4;
448   PetscScalar       dPD_dVr,dPD_dVi,dQD_dVr,dQD_dVi;
449   PetscScalar       dIDr_dVr,dIDr_dVi,dIDi_dVr,dIDi_dVi;
450 
451   PetscFunctionBegin;
452   PetscCall(MatZeroEntries(B));
453   PetscCall(DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet));
454   PetscCall(DMCompositeScatter(user->dmpgrid,X,Xgen,Xnet));
455 
456   PetscCall(VecGetArray(Xgen,&xgen));
457   PetscCall(VecGetArray(Xnet,&xnet));
458 
459   /* Generator subsystem */
460   for (i=0; i < ngen; i++) {
461     Eqp   = xgen[idx];
462     Edp   = xgen[idx+1];
463     delta = xgen[idx+2];
464     Id    = xgen[idx+4];
465     Iq    = xgen[idx+5];
466     Efd   = xgen[idx+6];
467 
468     /*    fgen[idx]   = (Eqp + (Xd[i] - Xdp[i])*Id - Efd)/Td0p[i]; */
469     row[0] = idx;
470     col[0] = idx;           col[1] = idx+4;          col[2] = idx+6;
471     val[0] = 1/ Td0p[i]; val[1] = (Xd[i] - Xdp[i])/ Td0p[i]; val[2] = -1/Td0p[i];
472 
473     PetscCall(MatSetValues(J,1,row,3,col,val,INSERT_VALUES));
474 
475     /*    fgen[idx+1] = (Edp - (Xq[i] - Xqp[i])*Iq)/Tq0p[i]; */
476     row[0] = idx + 1;
477     col[0] = idx + 1;       col[1] = idx+5;
478     val[0] = 1/Tq0p[i]; val[1] = -(Xq[i] - Xqp[i])/Tq0p[i];
479     PetscCall(MatSetValues(J,1,row,2,col,val,INSERT_VALUES));
480 
481     /*    fgen[idx+2] = - w + w_s; */
482     row[0] = idx + 2;
483     col[0] = idx + 2; col[1] = idx + 3;
484     val[0] = 0;       val[1] = -1;
485     PetscCall(MatSetValues(J,1,row,2,col,val,INSERT_VALUES));
486 
487     /*    fgen[idx+3] = (-TM[i] + Edp*Id + Eqp*Iq + (Xqp[i] - Xdp[i])*Id*Iq + D[i]*(w - w_s))/M[i]; */
488     row[0] = idx + 3;
489     col[0] = idx; col[1] = idx + 1; col[2] = idx + 3;       col[3] = idx + 4;                  col[4] = idx + 5;
490     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];
491     PetscCall(MatSetValues(J,1,row,5,col,val,INSERT_VALUES));
492 
493     Vr   = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */
494     Vi   = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */
495     PetscCall(ri2dq(Vr,Vi,delta,&Vd,&Vq));
496 
497     det = Rs[i]*Rs[i] + Xdp[i]*Xqp[i];
498 
499     Zdq_inv[0] = Rs[i]/det;
500     Zdq_inv[1] = Xqp[i]/det;
501     Zdq_inv[2] = -Xdp[i]/det;
502     Zdq_inv[3] = Rs[i]/det;
503 
504     dVd_dVr    = PetscSinScalar(delta); dVd_dVi = -PetscCosScalar(delta);
505     dVq_dVr    = PetscCosScalar(delta); dVq_dVi = PetscSinScalar(delta);
506     dVd_ddelta = Vr*PetscCosScalar(delta) + Vi*PetscSinScalar(delta);
507     dVq_ddelta = -Vr*PetscSinScalar(delta) + Vi*PetscCosScalar(delta);
508 
509     /*    fgen[idx+4] = Zdq_inv[0]*(-Edp + Vd) + Zdq_inv[1]*(-Eqp + Vq) + Id; */
510     row[0] = idx+4;
511     col[0] = idx;         col[1] = idx+1;        col[2] = idx + 2;
512     val[0] = -Zdq_inv[1]; val[1] = -Zdq_inv[0];  val[2] = Zdq_inv[0]*dVd_ddelta + Zdq_inv[1]*dVq_ddelta;
513     col[3] = idx + 4; col[4] = net_start+2*gbus[i];                     col[5] = net_start + 2*gbus[i]+1;
514     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;
515     PetscCall(MatSetValues(J,1,row,6,col,val,INSERT_VALUES));
516 
517     /*  fgen[idx+5] = Zdq_inv[2]*(-Edp + Vd) + Zdq_inv[3]*(-Eqp + Vq) + Iq; */
518     row[0] = idx+5;
519     col[0] = idx;         col[1] = idx+1;        col[2] = idx + 2;
520     val[0] = -Zdq_inv[3]; val[1] = -Zdq_inv[2];  val[2] = Zdq_inv[2]*dVd_ddelta + Zdq_inv[3]*dVq_ddelta;
521     col[3] = idx + 5; col[4] = net_start+2*gbus[i];                     col[5] = net_start + 2*gbus[i]+1;
522     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;
523     PetscCall(MatSetValues(J,1,row,6,col,val,INSERT_VALUES));
524 
525     dIGr_ddelta = Id*PetscCosScalar(delta) - Iq*PetscSinScalar(delta);
526     dIGi_ddelta = Id*PetscSinScalar(delta) + Iq*PetscCosScalar(delta);
527     dIGr_dId    = PetscSinScalar(delta);  dIGr_dIq = PetscCosScalar(delta);
528     dIGi_dId    = -PetscCosScalar(delta); dIGi_dIq = PetscSinScalar(delta);
529 
530     /* fnet[2*gbus[i]]   -= IGi; */
531     row[0] = net_start + 2*gbus[i];
532     col[0] = idx+2;        col[1] = idx + 4;   col[2] = idx + 5;
533     val[0] = -dIGi_ddelta; val[1] = -dIGi_dId; val[2] = -dIGi_dIq;
534     PetscCall(MatSetValues(J,1,row,3,col,val,INSERT_VALUES));
535 
536     /* fnet[2*gbus[i]+1]   -= IGr; */
537     row[0] = net_start + 2*gbus[i]+1;
538     col[0] = idx+2;        col[1] = idx + 4;   col[2] = idx + 5;
539     val[0] = -dIGr_ddelta; val[1] = -dIGr_dId; val[2] = -dIGr_dIq;
540     PetscCall(MatSetValues(J,1,row,3,col,val,INSERT_VALUES));
541 
542     Vm = PetscSqrtScalar(Vd*Vd + Vq*Vq);
543 
544     /*    fgen[idx+6] = (KE[i]*Efd + SE - VR)/TE[i]; */
545     /*    SE  = k1[i]*PetscExpScalar(k2[i]*Efd); */
546 
547     dSE_dEfd = k1[i]*k2[i]*PetscExpScalar(k2[i]*Efd);
548 
549     row[0] = idx + 6;
550     col[0] = idx + 6;                     col[1] = idx + 8;
551     val[0] = (KE[i] + dSE_dEfd)/TE[i];  val[1] = -1/TE[i];
552     PetscCall(MatSetValues(J,1,row,2,col,val,INSERT_VALUES));
553 
554     /* Exciter differential equations */
555 
556     /*    fgen[idx+7] = (RF - KF[i]*Efd/TF[i])/TF[i]; */
557     row[0] = idx + 7;
558     col[0] = idx + 6;       col[1] = idx + 7;
559     val[0] = (-KF[i]/TF[i])/TF[i];  val[1] = 1/TF[i];
560     PetscCall(MatSetValues(J,1,row,2,col,val,INSERT_VALUES));
561 
562     /*    fgen[idx+8] = (VR - KA[i]*RF + KA[i]*KF[i]*Efd/TF[i] - KA[i]*(Vref[i] - Vm))/TA[i]; */
563     /* Vm = (Vd^2 + Vq^2)^0.5; */
564     dVm_dVd    = Vd/Vm; dVm_dVq = Vq/Vm;
565     dVm_dVr    = dVm_dVd*dVd_dVr + dVm_dVq*dVq_dVr;
566     dVm_dVi    = dVm_dVd*dVd_dVi + dVm_dVq*dVq_dVi;
567     row[0]     = idx + 8;
568     col[0]     = idx + 6;           col[1] = idx + 7; col[2] = idx + 8;
569     val[0]     = (KA[i]*KF[i]/TF[i])/TA[i]; val[1] = -KA[i]/TA[i];  val[2] = 1/TA[i];
570     col[3]     = net_start + 2*gbus[i]; col[4] = net_start + 2*gbus[i]+1;
571     val[3]     = KA[i]*dVm_dVr/TA[i];         val[4] = KA[i]*dVm_dVi/TA[i];
572     PetscCall(MatSetValues(J,1,row,5,col,val,INSERT_VALUES));
573     idx        = idx + 9;
574   }
575 
576   for (i=0; i<nbus; i++) {
577     PetscCall(MatGetRow(user->Ybus,2*i,&ncols,&cols,&yvals));
578     row[0] = net_start + 2*i;
579     for (k=0; k<ncols; k++) {
580       col[k] = net_start + cols[k];
581       val[k] = yvals[k];
582     }
583     PetscCall(MatSetValues(J,1,row,ncols,col,val,INSERT_VALUES));
584     PetscCall(MatRestoreRow(user->Ybus,2*i,&ncols,&cols,&yvals));
585 
586     PetscCall(MatGetRow(user->Ybus,2*i+1,&ncols,&cols,&yvals));
587     row[0] = net_start + 2*i+1;
588     for (k=0; k<ncols; k++) {
589       col[k] = net_start + cols[k];
590       val[k] = yvals[k];
591     }
592     PetscCall(MatSetValues(J,1,row,ncols,col,val,INSERT_VALUES));
593     PetscCall(MatRestoreRow(user->Ybus,2*i+1,&ncols,&cols,&yvals));
594   }
595 
596   PetscCall(MatAssemblyBegin(J,MAT_FLUSH_ASSEMBLY));
597   PetscCall(MatAssemblyEnd(J,MAT_FLUSH_ASSEMBLY));
598 
599   PetscCall(VecGetArray(user->V0,&v0));
600   for (i=0; i < nload; i++) {
601     Vr      = xnet[2*lbus[i]]; /* Real part of load bus voltage */
602     Vi      = xnet[2*lbus[i]+1]; /* Imaginary part of the load bus voltage */
603     Vm      = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm; Vm4 = Vm2*Vm2;
604     Vm0     = PetscSqrtScalar(v0[2*lbus[i]]*v0[2*lbus[i]] + v0[2*lbus[i]+1]*v0[2*lbus[i]+1]);
605     PD      = QD = 0.0;
606     dPD_dVr = dPD_dVi = dQD_dVr = dQD_dVi = 0.0;
607     for (k=0; k < ld_nsegsp[i]; k++) {
608       PD      += ld_alphap[k]*PD0[i]*PetscPowScalar((Vm/Vm0),ld_betap[k]);
609       dPD_dVr += ld_alphap[k]*ld_betap[k]*PD0[i]*PetscPowScalar((1/Vm0),ld_betap[k])*Vr*PetscPowScalar(Vm,(ld_betap[k]-2));
610       dPD_dVi += ld_alphap[k]*ld_betap[k]*PD0[i]*PetscPowScalar((1/Vm0),ld_betap[k])*Vi*PetscPowScalar(Vm,(ld_betap[k]-2));
611     }
612     for (k=0; k < ld_nsegsq[i]; k++) {
613       QD      += ld_alphaq[k]*QD0[i]*PetscPowScalar((Vm/Vm0),ld_betaq[k]);
614       dQD_dVr += ld_alphaq[k]*ld_betaq[k]*QD0[i]*PetscPowScalar((1/Vm0),ld_betaq[k])*Vr*PetscPowScalar(Vm,(ld_betaq[k]-2));
615       dQD_dVi += ld_alphaq[k]*ld_betaq[k]*QD0[i]*PetscPowScalar((1/Vm0),ld_betaq[k])*Vi*PetscPowScalar(Vm,(ld_betaq[k]-2));
616     }
617 
618     /*    IDr = (PD*Vr + QD*Vi)/Vm2; */
619     /*    IDi = (-QD*Vr + PD*Vi)/Vm2; */
620 
621     dIDr_dVr = (dPD_dVr*Vr + dQD_dVr*Vi + PD)/Vm2 - ((PD*Vr + QD*Vi)*2*Vr)/Vm4;
622     dIDr_dVi = (dPD_dVi*Vr + dQD_dVi*Vi + QD)/Vm2 - ((PD*Vr + QD*Vi)*2*Vi)/Vm4;
623 
624     dIDi_dVr = (-dQD_dVr*Vr + dPD_dVr*Vi - QD)/Vm2 - ((-QD*Vr + PD*Vi)*2*Vr)/Vm4;
625     dIDi_dVi = (-dQD_dVi*Vr + dPD_dVi*Vi + PD)/Vm2 - ((-QD*Vr + PD*Vi)*2*Vi)/Vm4;
626 
627     /*    fnet[2*lbus[i]]   += IDi; */
628     row[0] = net_start + 2*lbus[i];
629     col[0] = net_start + 2*lbus[i];  col[1] = net_start + 2*lbus[i]+1;
630     val[0] = dIDi_dVr;               val[1] = dIDi_dVi;
631     PetscCall(MatSetValues(J,1,row,2,col,val,ADD_VALUES));
632     /*    fnet[2*lbus[i]+1] += IDr; */
633     row[0] = net_start + 2*lbus[i]+1;
634     col[0] = net_start + 2*lbus[i];  col[1] = net_start + 2*lbus[i]+1;
635     val[0] = dIDr_dVr;               val[1] = dIDr_dVi;
636     PetscCall(MatSetValues(J,1,row,2,col,val,ADD_VALUES));
637   }
638   PetscCall(VecRestoreArray(user->V0,&v0));
639 
640   PetscCall(VecRestoreArray(Xgen,&xgen));
641   PetscCall(VecRestoreArray(Xnet,&xnet));
642 
643   PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet));
644 
645   PetscCall(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY));
646   PetscCall(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY));
647   PetscFunctionReturn(0);
648 }
649 
650 /*
651    J = [I, 0
652         dg_dx, dg_dy]
653 */
654 PetscErrorCode AlgJacobian(SNES snes,Vec X,Mat A,Mat B,void *ctx)
655 {
656   Userctx        *user=(Userctx*)ctx;
657 
658   PetscFunctionBegin;
659   PetscCall(ResidualJacobian(snes,X,A,B,ctx));
660   PetscCall(MatSetOption(A,MAT_KEEP_NONZERO_PATTERN,PETSC_TRUE));
661   PetscCall(MatZeroRowsIS(A,user->is_diff,1.0,NULL,NULL));
662   PetscFunctionReturn(0);
663 }
664 
665 /*
666    J = [a*I-df_dx, -df_dy
667         dg_dx, dg_dy]
668 */
669 
670 PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,Userctx *user)
671 {
672   SNES           snes;
673   PetscScalar    atmp = (PetscScalar) a;
674   PetscInt       i,row;
675 
676   PetscFunctionBegin;
677   user->t = t;
678 
679   PetscCall(TSGetSNES(ts,&snes));
680   PetscCall(ResidualJacobian(snes,X,A,B,user));
681   for (i=0;i < ngen;i++) {
682     row = 9*i;
683     PetscCall(MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES));
684     row  = 9*i+1;
685     PetscCall(MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES));
686     row  = 9*i+2;
687     PetscCall(MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES));
688     row  = 9*i+3;
689     PetscCall(MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES));
690     row  = 9*i+6;
691     PetscCall(MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES));
692     row  = 9*i+7;
693     PetscCall(MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES));
694     row  = 9*i+8;
695     PetscCall(MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES));
696   }
697   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
698   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
699   PetscFunctionReturn(0);
700 }
701 
702 static PetscErrorCode CostIntegrand(TS ts,PetscReal t,Vec U,Vec R,Userctx *user)
703 {
704   PetscScalar       *r;
705   const PetscScalar *u;
706   PetscInt          idx;
707   Vec               Xgen,Xnet;
708   PetscScalar       *xgen;
709   PetscInt          i;
710 
711   PetscFunctionBegin;
712   PetscCall(DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet));
713   PetscCall(DMCompositeScatter(user->dmpgrid,U,Xgen,Xnet));
714 
715   PetscCall(VecGetArray(Xgen,&xgen));
716 
717   PetscCall(VecGetArrayRead(U,&u));
718   PetscCall(VecGetArray(R,&r));
719   r[0] = 0.;
720 
721   idx = 0;
722   for (i=0;i<ngen;i++) {
723     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);
724     idx  += 9;
725   }
726   PetscCall(VecRestoreArray(R,&r));
727   PetscCall(VecRestoreArrayRead(U,&u));
728   PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet));
729   PetscFunctionReturn(0);
730 }
731 
732 static PetscErrorCode MonitorUpdateQ(TS ts,PetscInt stepnum,PetscReal time,Vec X,void *ctx0)
733 {
734   Vec            C,*Y;
735   PetscInt       Nr;
736   PetscReal      h,theta;
737   Userctx        *ctx=(Userctx*)ctx0;
738 
739   PetscFunctionBegin;
740   theta = 0.5;
741   PetscCall(TSGetStages(ts,&Nr,&Y));
742   PetscCall(TSGetTimeStep(ts,&h));
743   PetscCall(VecDuplicate(ctx->vec_q,&C));
744   /* compute integrals */
745   if (stepnum>0) {
746     PetscCall(CostIntegrand(ts,time,X,C,ctx));
747     PetscCall(VecAXPY(ctx->vec_q,h*theta,C));
748     PetscCall(CostIntegrand(ts,time+h*theta,Y[0],C,ctx));
749     PetscCall(VecAXPY(ctx->vec_q,h*(1-theta),C));
750   }
751   PetscCall(VecDestroy(&C));
752   PetscFunctionReturn(0);
753 }
754 
755 int main(int argc,char **argv)
756 {
757   Userctx            user;
758   Vec                p;
759   PetscScalar        *x_ptr;
760   PetscMPIInt        size;
761   PetscInt           i;
762   KSP                ksp;
763   PC                 pc;
764   PetscInt           *idx2;
765   Tao                tao;
766   Vec                lowerb,upperb;
767 
768   PetscFunctionBeginUser;
769   PetscCall(PetscInitialize(&argc,&argv,"petscoptions",help));
770   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
771   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Only for sequential runs");
772 
773   PetscCall(VecCreateSeq(PETSC_COMM_WORLD,1,&user.vec_q));
774 
775   user.neqs_gen   = 9*ngen; /* # eqs. for generator subsystem */
776   user.neqs_net   = 2*nbus; /* # eqs. for network subsystem   */
777   user.neqs_pgrid = user.neqs_gen + user.neqs_net;
778 
779   /* Create indices for differential and algebraic equations */
780   PetscCall(PetscMalloc1(7*ngen,&idx2));
781   for (i=0; i<ngen; i++) {
782     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;
783     idx2[7*i+4] = 9*i+6; idx2[7*i+5] = 9*i+7; idx2[7*i+6] = 9*i+8;
784   }
785   PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,7*ngen,idx2,PETSC_COPY_VALUES,&user.is_diff));
786   PetscCall(ISComplement(user.is_diff,0,user.neqs_pgrid,&user.is_alg));
787   PetscCall(PetscFree(idx2));
788 
789   /* Set run time options */
790   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Transient stability fault options","");
791   {
792     user.tfaulton  = 1.0;
793     user.tfaultoff = 1.2;
794     user.Rfault    = 0.0001;
795     user.faultbus  = 8;
796     PetscCall(PetscOptionsReal("-tfaulton","","",user.tfaulton,&user.tfaulton,NULL));
797     PetscCall(PetscOptionsReal("-tfaultoff","","",user.tfaultoff,&user.tfaultoff,NULL));
798     PetscCall(PetscOptionsInt("-faultbus","","",user.faultbus,&user.faultbus,NULL));
799     user.t0        = 0.0;
800     user.tmax      = 1.5;
801     PetscCall(PetscOptionsReal("-t0","","",user.t0,&user.t0,NULL));
802     PetscCall(PetscOptionsReal("-tmax","","",user.tmax,&user.tmax,NULL));
803     user.freq_u    = 61.0;
804     user.freq_l    = 59.0;
805     user.pow       = 2;
806     PetscCall(PetscOptionsReal("-frequ","","",user.freq_u,&user.freq_u,NULL));
807     PetscCall(PetscOptionsReal("-freql","","",user.freq_l,&user.freq_l,NULL));
808     PetscCall(PetscOptionsInt("-pow","","",user.pow,&user.pow,NULL));
809 
810   }
811   PetscOptionsEnd();
812 
813   /* Create DMs for generator and network subsystems */
814   PetscCall(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,user.neqs_gen,1,1,NULL,&user.dmgen));
815   PetscCall(DMSetOptionsPrefix(user.dmgen,"dmgen_"));
816   PetscCall(DMSetFromOptions(user.dmgen));
817   PetscCall(DMSetUp(user.dmgen));
818   PetscCall(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,user.neqs_net,1,1,NULL,&user.dmnet));
819   PetscCall(DMSetOptionsPrefix(user.dmnet,"dmnet_"));
820   PetscCall(DMSetFromOptions(user.dmnet));
821   PetscCall(DMSetUp(user.dmnet));
822   /* Create a composite DM packer and add the two DMs */
823   PetscCall(DMCompositeCreate(PETSC_COMM_WORLD,&user.dmpgrid));
824   PetscCall(DMSetOptionsPrefix(user.dmpgrid,"pgrid_"));
825   PetscCall(DMCompositeAddDM(user.dmpgrid,user.dmgen));
826   PetscCall(DMCompositeAddDM(user.dmpgrid,user.dmnet));
827 
828   /* Create TAO solver and set desired solution method */
829   PetscCall(TaoCreate(PETSC_COMM_WORLD,&tao));
830   PetscCall(TaoSetType(tao,TAOBLMVM));
831   /*
832      Optimization starts
833   */
834   /* Set initial solution guess */
835   PetscCall(VecCreateSeq(PETSC_COMM_WORLD,3,&p));
836   PetscCall(VecGetArray(p,&x_ptr));
837   x_ptr[0] = PG[0]; x_ptr[1] = PG[1]; x_ptr[2] = PG[2];
838   PetscCall(VecRestoreArray(p,&x_ptr));
839 
840   PetscCall(TaoSetSolution(tao,p));
841   /* Set routine for function and gradient evaluation */
842   PetscCall(TaoSetObjective(tao,FormFunction,(void *)&user));
843   PetscCall(TaoSetGradient(tao,NULL,TaoDefaultComputeGradient,(void *)&user));
844 
845   /* Set bounds for the optimization */
846   PetscCall(VecDuplicate(p,&lowerb));
847   PetscCall(VecDuplicate(p,&upperb));
848   PetscCall(VecGetArray(lowerb,&x_ptr));
849   x_ptr[0] = 0.5; x_ptr[1] = 0.5; x_ptr[2] = 0.5;
850   PetscCall(VecRestoreArray(lowerb,&x_ptr));
851   PetscCall(VecGetArray(upperb,&x_ptr));
852   x_ptr[0] = 2.0; x_ptr[1] = 2.0; x_ptr[2] = 2.0;
853   PetscCall(VecRestoreArray(upperb,&x_ptr));
854   PetscCall(TaoSetVariableBounds(tao,lowerb,upperb));
855 
856   /* Check for any TAO command line options */
857   PetscCall(TaoSetFromOptions(tao));
858   PetscCall(TaoGetKSP(tao,&ksp));
859   if (ksp) {
860     PetscCall(KSPGetPC(ksp,&pc));
861     PetscCall(PCSetType(pc,PCNONE));
862   }
863 
864   /* SOLVE THE APPLICATION */
865   PetscCall(TaoSolve(tao));
866 
867   PetscCall(VecView(p,PETSC_VIEWER_STDOUT_WORLD));
868   /* Free TAO data structures */
869   PetscCall(TaoDestroy(&tao));
870   PetscCall(VecDestroy(&user.vec_q));
871   PetscCall(VecDestroy(&lowerb));
872   PetscCall(VecDestroy(&upperb));
873   PetscCall(VecDestroy(&p));
874   PetscCall(DMDestroy(&user.dmgen));
875   PetscCall(DMDestroy(&user.dmnet));
876   PetscCall(DMDestroy(&user.dmpgrid));
877   PetscCall(ISDestroy(&user.is_diff));
878   PetscCall(ISDestroy(&user.is_alg));
879   PetscCall(PetscFinalize());
880   return 0;
881 }
882 
883 /* ------------------------------------------------------------------ */
884 /*
885    FormFunction - Evaluates the function and corresponding gradient.
886 
887    Input Parameters:
888    tao - the Tao context
889    X   - the input vector
890    ptr - optional user-defined context, as set by TaoSetObjectiveAndGradient()
891 
892    Output Parameters:
893    f   - the newly evaluated function
894 */
895 PetscErrorCode FormFunction(Tao tao,Vec P,PetscReal *f,void *ctx0)
896 {
897   TS             ts;
898   SNES           snes_alg;
899   Userctx        *ctx = (Userctx*)ctx0;
900   Vec            X;
901   Mat            J;
902   /* sensitivity context */
903   PetscScalar    *x_ptr;
904   PetscViewer    Xview,Ybusview;
905   Vec            F_alg;
906   Vec            Xdot;
907   PetscInt       row_loc,col_loc;
908   PetscScalar    val;
909 
910   PetscCall(VecGetArrayRead(P,(const PetscScalar**)&x_ptr));
911   PG[0] = x_ptr[0];
912   PG[1] = x_ptr[1];
913   PG[2] = x_ptr[2];
914   PetscCall(VecRestoreArrayRead(P,(const PetscScalar**)&x_ptr));
915 
916   ctx->stepnum = 0;
917 
918   PetscCall(VecZeroEntries(ctx->vec_q));
919 
920   /* Read initial voltage vector and Ybus */
921   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"X.bin",FILE_MODE_READ,&Xview));
922   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"Ybus.bin",FILE_MODE_READ,&Ybusview));
923 
924   PetscCall(VecCreate(PETSC_COMM_WORLD,&ctx->V0));
925   PetscCall(VecSetSizes(ctx->V0,PETSC_DECIDE,ctx->neqs_net));
926   PetscCall(VecLoad(ctx->V0,Xview));
927 
928   PetscCall(MatCreate(PETSC_COMM_WORLD,&ctx->Ybus));
929   PetscCall(MatSetSizes(ctx->Ybus,PETSC_DECIDE,PETSC_DECIDE,ctx->neqs_net,ctx->neqs_net));
930   PetscCall(MatSetType(ctx->Ybus,MATBAIJ));
931   /*  PetscCall(MatSetBlockSize(ctx->Ybus,2)); */
932   PetscCall(MatLoad(ctx->Ybus,Ybusview));
933 
934   PetscCall(PetscViewerDestroy(&Xview));
935   PetscCall(PetscViewerDestroy(&Ybusview));
936 
937   PetscCall(DMCreateGlobalVector(ctx->dmpgrid,&X));
938 
939   PetscCall(MatCreate(PETSC_COMM_WORLD,&J));
940   PetscCall(MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,ctx->neqs_pgrid,ctx->neqs_pgrid));
941   PetscCall(MatSetFromOptions(J));
942   PetscCall(PreallocateJacobian(J,ctx));
943 
944   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
945      Create timestepping solver context
946      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
947   PetscCall(TSCreate(PETSC_COMM_WORLD,&ts));
948   PetscCall(TSSetProblemType(ts,TS_NONLINEAR));
949   PetscCall(TSSetType(ts,TSCN));
950   PetscCall(TSSetIFunction(ts,NULL,(TSIFunction) IFunction,ctx));
951   PetscCall(TSSetIJacobian(ts,J,J,(TSIJacobian)IJacobian,ctx));
952   PetscCall(TSSetApplicationContext(ts,ctx));
953 
954   PetscCall(TSMonitorSet(ts,MonitorUpdateQ,ctx,NULL));
955   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
956      Set initial conditions
957    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
958   PetscCall(SetInitialGuess(X,ctx));
959 
960   PetscCall(VecDuplicate(X,&F_alg));
961   PetscCall(SNESCreate(PETSC_COMM_WORLD,&snes_alg));
962   PetscCall(SNESSetFunction(snes_alg,F_alg,AlgFunction,ctx));
963   PetscCall(MatZeroEntries(J));
964   PetscCall(SNESSetJacobian(snes_alg,J,J,AlgJacobian,ctx));
965   PetscCall(SNESSetOptionsPrefix(snes_alg,"alg_"));
966   PetscCall(SNESSetFromOptions(snes_alg));
967   ctx->alg_flg = PETSC_TRUE;
968   /* Solve the algebraic equations */
969   PetscCall(SNESSolve(snes_alg,NULL,X));
970 
971   /* Just to set up the Jacobian structure */
972   PetscCall(VecDuplicate(X,&Xdot));
973   PetscCall(IJacobian(ts,0.0,X,Xdot,0.0,J,J,ctx));
974   PetscCall(VecDestroy(&Xdot));
975 
976   ctx->stepnum++;
977 
978   PetscCall(TSSetTimeStep(ts,0.01));
979   PetscCall(TSSetMaxTime(ts,ctx->tmax));
980   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP));
981   PetscCall(TSSetFromOptions(ts));
982 
983   /* Prefault period */
984   ctx->alg_flg = PETSC_FALSE;
985   PetscCall(TSSetTime(ts,0.0));
986   PetscCall(TSSetMaxTime(ts,ctx->tfaulton));
987   PetscCall(TSSolve(ts,X));
988 
989   /* Create the nonlinear solver for solving the algebraic system */
990   /* Note that although the algebraic system needs to be solved only for
991      Idq and V, we reuse the entire system including xgen. The xgen
992      variables are held constant by setting their residuals to 0 and
993      putting a 1 on the Jacobian diagonal for xgen rows
994   */
995   PetscCall(MatZeroEntries(J));
996 
997   /* Apply disturbance - resistive fault at ctx->faultbus */
998   /* This is done by adding shunt conductance to the diagonal location
999      in the Ybus matrix */
1000   row_loc = 2*ctx->faultbus; col_loc = 2*ctx->faultbus+1; /* Location for G */
1001   val     = 1/ctx->Rfault;
1002   PetscCall(MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES));
1003   row_loc = 2*ctx->faultbus+1; col_loc = 2*ctx->faultbus; /* Location for G */
1004   val     = 1/ctx->Rfault;
1005   PetscCall(MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES));
1006 
1007   PetscCall(MatAssemblyBegin(ctx->Ybus,MAT_FINAL_ASSEMBLY));
1008   PetscCall(MatAssemblyEnd(ctx->Ybus,MAT_FINAL_ASSEMBLY));
1009 
1010   ctx->alg_flg = PETSC_TRUE;
1011   /* Solve the algebraic equations */
1012   PetscCall(SNESSolve(snes_alg,NULL,X));
1013 
1014   ctx->stepnum++;
1015 
1016   /* Disturbance period */
1017   ctx->alg_flg = PETSC_FALSE;
1018   PetscCall(TSSetTime(ts,ctx->tfaulton));
1019   PetscCall(TSSetMaxTime(ts,ctx->tfaultoff));
1020   PetscCall(TSSolve(ts,X));
1021 
1022   /* Remove the fault */
1023   row_loc = 2*ctx->faultbus; col_loc = 2*ctx->faultbus+1;
1024   val     = -1/ctx->Rfault;
1025   PetscCall(MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES));
1026   row_loc = 2*ctx->faultbus+1; col_loc = 2*ctx->faultbus;
1027   val     = -1/ctx->Rfault;
1028   PetscCall(MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES));
1029 
1030   PetscCall(MatAssemblyBegin(ctx->Ybus,MAT_FINAL_ASSEMBLY));
1031   PetscCall(MatAssemblyEnd(ctx->Ybus,MAT_FINAL_ASSEMBLY));
1032 
1033   PetscCall(MatZeroEntries(J));
1034 
1035   ctx->alg_flg = PETSC_TRUE;
1036 
1037   /* Solve the algebraic equations */
1038   PetscCall(SNESSolve(snes_alg,NULL,X));
1039 
1040   ctx->stepnum++;
1041 
1042   /* Post-disturbance period */
1043   ctx->alg_flg = PETSC_TRUE;
1044   PetscCall(TSSetTime(ts,ctx->tfaultoff));
1045   PetscCall(TSSetMaxTime(ts,ctx->tmax));
1046   PetscCall(TSSolve(ts,X));
1047 
1048   PetscCall(VecGetArray(ctx->vec_q,&x_ptr));
1049   *f   = x_ptr[0];
1050   PetscCall(VecRestoreArray(ctx->vec_q,&x_ptr));
1051 
1052   PetscCall(MatDestroy(&ctx->Ybus));
1053   PetscCall(VecDestroy(&ctx->V0));
1054   PetscCall(SNESDestroy(&snes_alg));
1055   PetscCall(VecDestroy(&F_alg));
1056   PetscCall(MatDestroy(&J));
1057   PetscCall(VecDestroy(&X));
1058   PetscCall(TSDestroy(&ts));
1059 
1060   return 0;
1061 }
1062 
1063 /*TEST
1064 
1065   build:
1066       requires: double !complex !defined(USE_64BIT_INDICES)
1067 
1068    test:
1069       args: -viewer_binary_skip_info -tao_monitor -tao_gttol .2
1070       localrunfiles: petscoptions X.bin Ybus.bin
1071 
1072 TEST*/
1073