xref: /petsc/src/ts/tutorials/power_grid/stability_9bus/ex9bus.c (revision 6a98f8dc3f2c9149905a87dc2e9d0fedaf64e09a)
1 
2 static char help[] = "Power grid stability analysis of WECC 9 bus system.\n\
3 This example is based on the 9-bus (node) example given in the book Power\n\
4 Systems Dynamics and Stability (Chapter 7) by P. Sauer and M. A. Pai.\n\
5 The power grid in this example consists of 9 buses (nodes), 3 generators,\n\
6 3 loads, and 9 transmission lines. The network equations are written\n\
7 in current balance form using rectangular coordiantes.\n\n";
8 
9 /*
10    The equations for the stability analysis are described by the DAE
11 
12    \dot{x} = f(x,y,t)
13      0     = g(x,y,t)
14 
15    where the generators are described by differential equations, while the algebraic
16    constraints define the network equations.
17 
18    The generators are modeled with a 4th order differential equation describing the electrical
19    and mechanical dynamics. Each generator also has an exciter system modeled by 3rd order
20    diff. eqns. describing the exciter, voltage regulator, and the feedback stabilizer
21    mechanism.
22 
23    The network equations are described by nodal current balance equations.
24     I(x,y) - Y*V = 0
25 
26    where:
27     I(x,y) is the current injected from generators and loads.
28       Y    is the admittance matrix, and
29       V    is the voltage vector
30 */
31 
32 /*
33    Include "petscts.h" so that we can use TS solvers.  Note that this
34    file automatically includes:
35      petscsys.h       - base PETSc routines   petscvec.h - vectors
36      petscmat.h - matrices
37      petscis.h     - index sets            petscksp.h - Krylov subspace methods
38      petscviewer.h - viewers               petscpc.h  - preconditioners
39      petscksp.h   - linear solvers
40 */
41 
42 #include <petscts.h>
43 #include <petscdm.h>
44 #include <petscdmda.h>
45 #include <petscdmcomposite.h>
46 
47 #define freq 60
48 #define w_s (2*PETSC_PI*freq)
49 
50 /* Sizes and indices */
51 const PetscInt nbus    = 9; /* Number of network buses */
52 const PetscInt ngen    = 3; /* Number of generators */
53 const PetscInt nload   = 3; /* Number of loads */
54 const PetscInt gbus[3] = {0,1,2}; /* Buses at which generators are incident */
55 const PetscInt lbus[3] = {4,5,7}; /* Buses at which loads are incident */
56 
57 /* Generator real and reactive powers (found via loadflow) */
58 const PetscScalar PG[3] = {0.716786142395021,1.630000000000000,0.850000000000000};
59 const PetscScalar QG[3] = {0.270702180178785,0.066120127797275,-0.108402221791588};
60 /* Generator constants */
61 const PetscScalar H[3]    = {23.64,6.4,3.01};   /* Inertia constant */
62 const PetscScalar Rs[3]   = {0.0,0.0,0.0}; /* Stator Resistance */
63 const PetscScalar Xd[3]   = {0.146,0.8958,1.3125};  /* d-axis reactance */
64 const PetscScalar Xdp[3]  = {0.0608,0.1198,0.1813}; /* d-axis transient reactance */
65 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 */
66 const PetscScalar Xqp[3]  = {0.0969,0.1969,0.25}; /* q-axis transient reactance */
67 const PetscScalar Td0p[3] = {8.96,6.0,5.89}; /* d-axis open circuit time constant */
68 const PetscScalar Tq0p[3] = {0.31,0.535,0.6}; /* q-axis open circuit time constant */
69 PetscScalar M[3]; /* M = 2*H/w_s */
70 PetscScalar D[3]; /* D = 0.1*M */
71 
72 PetscScalar TM[3]; /* Mechanical Torque */
73 /* Exciter system constants */
74 const PetscScalar KA[3] = {20.0,20.0,20.0};  /* Voltage regulartor gain constant */
75 const PetscScalar TA[3] = {0.2,0.2,0.2};     /* Voltage regulator time constant */
76 const PetscScalar KE[3] = {1.0,1.0,1.0};     /* Exciter gain constant */
77 const PetscScalar TE[3] = {0.314,0.314,0.314}; /* Exciter time constant */
78 const PetscScalar KF[3] = {0.063,0.063,0.063};  /* Feedback stabilizer gain constant */
79 const PetscScalar TF[3] = {0.35,0.35,0.35};    /* Feedback stabilizer time constant */
80 const PetscScalar k1[3] = {0.0039,0.0039,0.0039};
81 const PetscScalar k2[3] = {1.555,1.555,1.555};  /* k1 and k2 for calculating the saturation function SE = k1*exp(k2*Efd) */
82 const PetscScalar VRMIN[3] = {-4.0,-4.0,-4.0};
83 const PetscScalar VRMAX[3] = {7.0,7.0,7.0};
84 PetscInt VRatmin[3];
85 PetscInt VRatmax[3];
86 
87 PetscScalar Vref[3];
88 /* Load constants
89   We use a composite load model that describes the load and reactive powers at each time instant as follows
90   P(t) = \sum\limits_{i=0}^ld_nsegsp \ld_alphap_i*P_D0(\frac{V_m(t)}{V_m0})^\ld_betap_i
91   Q(t) = \sum\limits_{i=0}^ld_nsegsq \ld_alphaq_i*Q_D0(\frac{V_m(t)}{V_m0})^\ld_betaq_i
92   where
93     ld_nsegsp,ld_nsegsq - Number of individual load models for real and reactive power loads
94     ld_alphap,ld_alphap - Percentage contribution (weights) or loads
95     P_D0                - Real power load
96     Q_D0                - Reactive power load
97     V_m(t)              - Voltage magnitude at time t
98     V_m0                - Voltage magnitude at t = 0
99     ld_betap, ld_betaq  - exponents describing the load model for real and reactive part
100 
101     Note: All loads have the same characteristic currently.
102 */
103 const PetscScalar PD0[3] = {1.25,0.9,1.0};
104 const PetscScalar QD0[3] = {0.5,0.3,0.35};
105 const PetscInt    ld_nsegsp[3] = {3,3,3};
106 const PetscScalar ld_alphap[3] = {1.0,0.0,0.0};
107 const PetscScalar ld_betap[3]  = {2.0,1.0,0.0};
108 const PetscInt    ld_nsegsq[3] = {3,3,3};
109 const PetscScalar ld_alphaq[3] = {1.0,0.0,0.0};
110 const PetscScalar ld_betaq[3]  = {2.0,1.0,0.0};
111 
112 typedef struct {
113   DM          dmgen, dmnet; /* DMs to manage generator and network subsystem */
114   DM          dmpgrid; /* Composite DM to manage the entire power grid */
115   Mat         Ybus; /* Network admittance matrix */
116   Vec         V0;  /* Initial voltage vector (Power flow solution) */
117   PetscReal   tfaulton,tfaultoff; /* Fault on and off times */
118   PetscInt    faultbus; /* Fault bus */
119   PetscScalar Rfault;
120   PetscReal   t0,tmax;
121   PetscInt    neqs_gen,neqs_net,neqs_pgrid;
122   Mat         Sol; /* Matrix to save solution at each time step */
123   PetscInt    stepnum;
124   PetscReal   t;
125   SNES        snes_alg;
126   IS          is_diff; /* indices for differential equations */
127   IS          is_alg; /* indices for algebraic equations */
128   PetscBool   setisdiff; /* TS computes truncation error based only on the differential variables */
129   PetscBool   semiexplicit; /* If the flag is set then a semi-explicit method is used using TSRK */
130 } Userctx;
131 
132 /*
133   The first two events are for fault on and off, respectively. The following events are
134   to check the min/max limits on the state variable VR. A non windup limiter is used for
135   the VR limits.
136 */
137 PetscErrorCode EventFunction(TS ts,PetscReal t,Vec X,PetscScalar *fvalue,void *ctx)
138 {
139   Userctx        *user=(Userctx*)ctx;
140   Vec            Xgen,Xnet;
141   PetscInt       i,idx=0;
142   const PetscScalar *xgen,*xnet;
143   PetscErrorCode ierr;
144   PetscScalar    Efd,RF,VR,Vr,Vi,Vm;
145 
146   PetscFunctionBegin;
147 
148   ierr = DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);CHKERRQ(ierr);
149   ierr = DMCompositeScatter(user->dmpgrid,X,Xgen,Xnet);CHKERRQ(ierr);
150 
151   ierr = VecGetArrayRead(Xgen,&xgen);CHKERRQ(ierr);
152   ierr = VecGetArrayRead(Xnet,&xnet);CHKERRQ(ierr);
153 
154   /* Event for fault-on time */
155   fvalue[0] = t - user->tfaulton;
156   /* Event for fault-off time */
157   fvalue[1] = t - user->tfaultoff;
158 
159   for (i=0; i < ngen; i++) {
160     Efd   = xgen[idx+6];
161     RF    = xgen[idx+7];
162     VR    = xgen[idx+8];
163 
164     Vr = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */
165     Vi = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */
166     Vm = PetscSqrtScalar(Vr*Vr + Vi*Vi);
167 
168     if (!VRatmax[i]) {
169       fvalue[2+2*i] = VRMAX[i] - VR;
170     } else {
171       fvalue[2+2*i] = (VR - KA[i]*RF + KA[i]*KF[i]*Efd/TF[i] - KA[i]*(Vref[i] - Vm))/TA[i];
172     }
173     if (!VRatmin[i]) {
174       fvalue[2+2*i+1] = VRMIN[i] - VR;
175     } else {
176       fvalue[2+2*i+1] = (VR - KA[i]*RF + KA[i]*KF[i]*Efd/TF[i] - KA[i]*(Vref[i] - Vm))/TA[i];
177     }
178     idx = idx+9;
179   }
180   ierr = VecRestoreArrayRead(Xgen,&xgen);CHKERRQ(ierr);
181   ierr = VecRestoreArrayRead(Xnet,&xnet);CHKERRQ(ierr);
182 
183   ierr = DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);CHKERRQ(ierr);
184 
185   PetscFunctionReturn(0);
186 }
187 
188 PetscErrorCode PostEventFunction(TS ts,PetscInt nevents,PetscInt event_list[],PetscReal t,Vec X,PetscBool forwardsolve,void* ctx)
189 {
190   Userctx *user=(Userctx*)ctx;
191   Vec      Xgen,Xnet;
192   PetscScalar *xgen,*xnet;
193   PetscInt row_loc,col_loc;
194   PetscScalar val;
195   PetscErrorCode ierr;
196   PetscInt i,idx=0,event_num;
197   PetscScalar fvalue;
198   PetscScalar Efd, RF, VR;
199   PetscScalar Vr,Vi,Vm;
200 
201   PetscFunctionBegin;
202 
203   ierr = DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);CHKERRQ(ierr);
204   ierr = DMCompositeScatter(user->dmpgrid,X,Xgen,Xnet);CHKERRQ(ierr);
205 
206   ierr = VecGetArray(Xgen,&xgen);CHKERRQ(ierr);
207   ierr = VecGetArray(Xnet,&xnet);CHKERRQ(ierr);
208 
209   for (i=0; i < nevents; i++) {
210     if (event_list[i] == 0) {
211       /* Apply disturbance - resistive fault at user->faultbus */
212       /* This is done by adding shunt conductance to the diagonal location
213 	 in the Ybus matrix */
214       row_loc = 2*user->faultbus; col_loc = 2*user->faultbus+1; /* Location for G */
215       val     = 1/user->Rfault;
216       ierr    = MatSetValues(user->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);CHKERRQ(ierr);
217       row_loc = 2*user->faultbus+1; col_loc = 2*user->faultbus; /* Location for G */
218       val     = 1/user->Rfault;
219       ierr    = MatSetValues(user->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);CHKERRQ(ierr);
220 
221       ierr = MatAssemblyBegin(user->Ybus,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
222       ierr = MatAssemblyEnd(user->Ybus,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
223 
224       /* Solve the algebraic equations */
225       ierr = SNESSolve(user->snes_alg,NULL,X);CHKERRQ(ierr);
226     } else if(event_list[i] == 1) {
227       /* Remove the fault */
228       row_loc = 2*user->faultbus; col_loc = 2*user->faultbus+1;
229       val     = -1/user->Rfault;
230       ierr    = MatSetValues(user->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);CHKERRQ(ierr);
231       row_loc = 2*user->faultbus+1; col_loc = 2*user->faultbus;
232       val     = -1/user->Rfault;
233       ierr    = MatSetValues(user->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);CHKERRQ(ierr);
234 
235       ierr = MatAssemblyBegin(user->Ybus,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
236       ierr = MatAssemblyEnd(user->Ybus,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
237 
238       /* Solve the algebraic equations */
239       ierr = SNESSolve(user->snes_alg,NULL,X);CHKERRQ(ierr);
240 
241       /* Check the VR derivatives and reset flags if needed */
242       for (i=0; i < ngen; i++) {
243 	Efd   = xgen[idx+6];
244 	RF    = xgen[idx+7];
245 	VR    = xgen[idx+8];
246 
247 	Vr = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */
248 	Vi = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */
249 	Vm = PetscSqrtScalar(Vr*Vr + Vi*Vi);
250 
251 	if (VRatmax[i]) {
252 	  fvalue = (VR - KA[i]*RF + KA[i]*KF[i]*Efd/TF[i] - KA[i]*(Vref[i] - Vm))/TA[i];
253 	  if (fvalue < 0) {
254 	    VRatmax[i] = 0;
255 	    ierr = PetscPrintf(PETSC_COMM_SELF,"VR[%d]: dVR_dt went negative on fault clearing at time %g\n",i,t);CHKERRQ(ierr);
256 	  }
257 	}
258 	if (VRatmin[i]) {
259 	  fvalue =  (VR - KA[i]*RF + KA[i]*KF[i]*Efd/TF[i] - KA[i]*(Vref[i] - Vm))/TA[i];
260 
261 	  if(fvalue > 0) {
262 	    VRatmin[i] = 0;
263 	    ierr = PetscPrintf(PETSC_COMM_SELF,"VR[%d]: dVR_dt went positive on fault clearing at time %g\n",i,t);CHKERRQ(ierr);
264 	  }
265 	}
266 	idx = idx+9;
267       }
268     } else {
269       idx = (event_list[i]-2)/2;
270       event_num = (event_list[i]-2)%2;
271       if (event_num == 0) { /* Max VR */
272 	if (!VRatmax[idx]) {
273 	  VRatmax[idx] = 1;
274 	  ierr = PetscPrintf(PETSC_COMM_SELF,"VR[%d]: hit upper limit at time %g\n",idx,t);CHKERRQ(ierr);
275 	}
276 	else {
277 	  VRatmax[idx] = 0;
278 	  ierr = PetscPrintf(PETSC_COMM_SELF,"VR[%d]: freeing variable as dVR_dt is negative at time %g\n",idx,t);CHKERRQ(ierr);
279 	}
280       } else {
281 	if (!VRatmin[idx]) {
282 	  VRatmin[idx] = 1;
283 	  ierr = PetscPrintf(PETSC_COMM_SELF,"VR[%d]: hit lower limit at time %g\n",idx,t);CHKERRQ(ierr);
284 	}
285 	else {
286 	  VRatmin[idx] = 0;
287 	  ierr = PetscPrintf(PETSC_COMM_SELF,"VR[%d]: freeing variable as dVR_dt is positive at time %g\n",idx,t);CHKERRQ(ierr);
288 	}
289       }
290     }
291   }
292 
293   ierr = VecRestoreArray(Xgen,&xgen);CHKERRQ(ierr);
294   ierr = VecRestoreArray(Xnet,&xnet);CHKERRQ(ierr);
295 
296   ierr = DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);CHKERRQ(ierr);
297 
298   PetscFunctionReturn(0);
299 }
300 
301 /* Converts from machine frame (dq) to network (phase a real,imag) reference frame */
302 PetscErrorCode dq2ri(PetscScalar Fd,PetscScalar Fq,PetscScalar delta,PetscScalar *Fr, PetscScalar *Fi)
303 {
304   PetscFunctionBegin;
305   *Fr =  Fd*PetscSinScalar(delta) + Fq*PetscCosScalar(delta);
306   *Fi = -Fd*PetscCosScalar(delta) + Fq*PetscSinScalar(delta);
307   PetscFunctionReturn(0);
308 }
309 
310 /* Converts from network frame ([phase a real,imag) to machine (dq) reference frame */
311 PetscErrorCode ri2dq(PetscScalar Fr,PetscScalar Fi,PetscScalar delta,PetscScalar *Fd, PetscScalar *Fq)
312 {
313   PetscFunctionBegin;
314   *Fd =  Fr*PetscSinScalar(delta) - Fi*PetscCosScalar(delta);
315   *Fq =  Fr*PetscCosScalar(delta) + Fi*PetscSinScalar(delta);
316   PetscFunctionReturn(0);
317 }
318 
319 /* Saves the solution at each time to a matrix */
320 PetscErrorCode SaveSolution(TS ts)
321 {
322   PetscErrorCode    ierr;
323   Userctx           *user;
324   Vec               X;
325   const PetscScalar *x;
326   PetscScalar       *mat;
327   PetscInt          idx;
328   PetscReal         t;
329 
330   PetscFunctionBegin;
331   ierr     = TSGetApplicationContext(ts,&user);CHKERRQ(ierr);
332   ierr     = TSGetTime(ts,&t);CHKERRQ(ierr);
333   ierr     = TSGetSolution(ts,&X);CHKERRQ(ierr);
334   idx      = user->stepnum*(user->neqs_pgrid+1);
335   ierr     = MatDenseGetArray(user->Sol,&mat);CHKERRQ(ierr);
336   ierr     = VecGetArrayRead(X,&x);CHKERRQ(ierr);
337   mat[idx] = t;
338   ierr     = PetscArraycpy(mat+idx+1,x,user->neqs_pgrid);CHKERRQ(ierr);
339   ierr     = MatDenseRestoreArray(user->Sol,&mat);CHKERRQ(ierr);
340   ierr     = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
341   user->stepnum++;
342   PetscFunctionReturn(0);
343 }
344 
345 PetscErrorCode SetInitialGuess(Vec X,Userctx *user)
346 {
347   PetscErrorCode ierr;
348   Vec            Xgen,Xnet;
349   PetscScalar    *xgen,*xnet;
350   PetscInt       i,idx=0;
351   PetscScalar    Vr,Vi,IGr,IGi,Vm,Vm2;
352   PetscScalar    Eqp,Edp,delta;
353   PetscScalar    Efd,RF,VR; /* Exciter variables */
354   PetscScalar    Id,Iq;  /* Generator dq axis currents */
355   PetscScalar    theta,Vd,Vq,SE;
356 
357   PetscFunctionBegin;
358   M[0] = 2*H[0]/w_s; M[1] = 2*H[1]/w_s; M[2] = 2*H[2]/w_s;
359   D[0] = 0.1*M[0]; D[1] = 0.1*M[1]; D[2] = 0.1*M[2];
360 
361   ierr = DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);CHKERRQ(ierr);
362 
363   /* Network subsystem initialization */
364   ierr = VecCopy(user->V0,Xnet);CHKERRQ(ierr);
365 
366   /* Generator subsystem initialization */
367   ierr = VecGetArray(Xgen,&xgen);CHKERRQ(ierr);
368   ierr = VecGetArray(Xnet,&xnet);CHKERRQ(ierr);
369 
370   for (i=0; i < ngen; i++) {
371     Vr  = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */
372     Vi  = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */
373     Vm  = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm;
374     IGr = (Vr*PG[i] + Vi*QG[i])/Vm2;
375     IGi = (Vi*PG[i] - Vr*QG[i])/Vm2;
376 
377     delta = PetscAtan2Real(Vi+Xq[i]*IGr,Vr-Xq[i]*IGi); /* Machine angle */
378 
379     theta = PETSC_PI/2.0 - delta;
380 
381     Id = IGr*PetscCosScalar(theta) - IGi*PetscSinScalar(theta); /* d-axis stator current */
382     Iq = IGr*PetscSinScalar(theta) + IGi*PetscCosScalar(theta); /* q-axis stator current */
383 
384     Vd = Vr*PetscCosScalar(theta) - Vi*PetscSinScalar(theta);
385     Vq = Vr*PetscSinScalar(theta) + Vi*PetscCosScalar(theta);
386 
387     Edp = Vd + Rs[i]*Id - Xqp[i]*Iq; /* d-axis transient EMF */
388     Eqp = Vq + Rs[i]*Iq + Xdp[i]*Id; /* q-axis transient EMF */
389 
390     TM[i] = PG[i];
391 
392     /* The generator variables are ordered as [Eqp,Edp,delta,w,Id,Iq] */
393     xgen[idx]   = Eqp;
394     xgen[idx+1] = Edp;
395     xgen[idx+2] = delta;
396     xgen[idx+3] = w_s;
397 
398     idx = idx + 4;
399 
400     xgen[idx]   = Id;
401     xgen[idx+1] = Iq;
402 
403     idx = idx + 2;
404 
405     /* Exciter */
406     Efd = Eqp + (Xd[i] - Xdp[i])*Id;
407     SE  = k1[i]*PetscExpScalar(k2[i]*Efd);
408     VR  =  KE[i]*Efd + SE;
409     RF  =  KF[i]*Efd/TF[i];
410 
411     xgen[idx]   = Efd;
412     xgen[idx+1] = RF;
413     xgen[idx+2] = VR;
414 
415     Vref[i] = Vm + (VR/KA[i]);
416 
417     VRatmin[i] = VRatmax[i] = 0;
418 
419     idx = idx + 3;
420   }
421 
422   ierr = VecRestoreArray(Xgen,&xgen);CHKERRQ(ierr);
423   ierr = VecRestoreArray(Xnet,&xnet);CHKERRQ(ierr);
424 
425   /* ierr = VecView(Xgen,0);CHKERRQ(ierr); */
426   ierr = DMCompositeGather(user->dmpgrid,INSERT_VALUES,X,Xgen,Xnet);CHKERRQ(ierr);
427   ierr = DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);CHKERRQ(ierr);
428   PetscFunctionReturn(0);
429 }
430 
431 /* Computes F = [f(x,y);g(x,y)] */
432 PetscErrorCode ResidualFunction(Vec X, Vec F, Userctx *user)
433 {
434   PetscErrorCode ierr;
435   Vec            Xgen,Xnet,Fgen,Fnet;
436   PetscScalar    *xgen,*xnet,*fgen,*fnet;
437   PetscInt       i,idx=0;
438   PetscScalar    Vr,Vi,Vm,Vm2;
439   PetscScalar    Eqp,Edp,delta,w; /* Generator variables */
440   PetscScalar    Efd,RF,VR; /* Exciter variables */
441   PetscScalar    Id,Iq;  /* Generator dq axis currents */
442   PetscScalar    Vd,Vq,SE;
443   PetscScalar    IGr,IGi,IDr,IDi;
444   PetscScalar    Zdq_inv[4],det;
445   PetscScalar    PD,QD,Vm0,*v0;
446   PetscInt       k;
447 
448   PetscFunctionBegin;
449   ierr = VecZeroEntries(F);CHKERRQ(ierr);
450   ierr = DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);CHKERRQ(ierr);
451   ierr = DMCompositeGetLocalVectors(user->dmpgrid,&Fgen,&Fnet);CHKERRQ(ierr);
452   ierr = DMCompositeScatter(user->dmpgrid,X,Xgen,Xnet);CHKERRQ(ierr);
453   ierr = DMCompositeScatter(user->dmpgrid,F,Fgen,Fnet);CHKERRQ(ierr);
454 
455   /* Network current balance residual IG + Y*V + IL = 0. Only YV is added here.
456      The generator current injection, IG, and load current injection, ID are added later
457   */
458   /* Note that the values in Ybus are stored assuming the imaginary current balance
459      equation is ordered first followed by real current balance equation for each bus.
460      Thus imaginary current contribution goes in location 2*i, and
461      real current contribution in 2*i+1
462   */
463   ierr = MatMult(user->Ybus,Xnet,Fnet);CHKERRQ(ierr);
464 
465   ierr = VecGetArray(Xgen,&xgen);CHKERRQ(ierr);
466   ierr = VecGetArray(Xnet,&xnet);CHKERRQ(ierr);
467   ierr = VecGetArray(Fgen,&fgen);CHKERRQ(ierr);
468   ierr = VecGetArray(Fnet,&fnet);CHKERRQ(ierr);
469 
470   /* Generator subsystem */
471   for (i=0; i < ngen; i++) {
472     Eqp   = xgen[idx];
473     Edp   = xgen[idx+1];
474     delta = xgen[idx+2];
475     w     = xgen[idx+3];
476     Id    = xgen[idx+4];
477     Iq    = xgen[idx+5];
478     Efd   = xgen[idx+6];
479     RF    = xgen[idx+7];
480     VR    = xgen[idx+8];
481 
482     /* Generator differential equations */
483     fgen[idx]   = (-Eqp - (Xd[i] - Xdp[i])*Id + Efd)/Td0p[i];
484     fgen[idx+1] = (-Edp + (Xq[i] - Xqp[i])*Iq)/Tq0p[i];
485     fgen[idx+2] = w - w_s;
486     fgen[idx+3] = (TM[i] - Edp*Id - Eqp*Iq - (Xqp[i] - Xdp[i])*Id*Iq - D[i]*(w - w_s))/M[i];
487 
488     Vr = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */
489     Vi = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */
490 
491     ierr = ri2dq(Vr,Vi,delta,&Vd,&Vq);CHKERRQ(ierr);
492     /* Algebraic equations for stator currents */
493     det = Rs[i]*Rs[i] + Xdp[i]*Xqp[i];
494 
495     Zdq_inv[0] = Rs[i]/det;
496     Zdq_inv[1] = Xqp[i]/det;
497     Zdq_inv[2] = -Xdp[i]/det;
498     Zdq_inv[3] = Rs[i]/det;
499 
500     fgen[idx+4] = Zdq_inv[0]*(-Edp + Vd) + Zdq_inv[1]*(-Eqp + Vq) + Id;
501     fgen[idx+5] = Zdq_inv[2]*(-Edp + Vd) + Zdq_inv[3]*(-Eqp + Vq) + Iq;
502 
503     /* Add generator current injection to network */
504     ierr = dq2ri(Id,Iq,delta,&IGr,&IGi);CHKERRQ(ierr);
505 
506     fnet[2*gbus[i]]   -= IGi;
507     fnet[2*gbus[i]+1] -= IGr;
508 
509     Vm = PetscSqrtScalar(Vd*Vd + Vq*Vq);
510 
511     SE = k1[i]*PetscExpScalar(k2[i]*Efd);
512 
513     /* Exciter differential equations */
514     fgen[idx+6] = (-KE[i]*Efd - SE + VR)/TE[i];
515     fgen[idx+7] = (-RF + KF[i]*Efd/TF[i])/TF[i];
516     if(VRatmax[i]) fgen[idx+8] = VR - VRMAX[i];
517     else if(VRatmin[i]) fgen[idx+8] = VRMIN[i] - VR;
518     else fgen[idx+8] = (-VR + KA[i]*RF - KA[i]*KF[i]*Efd/TF[i] + KA[i]*(Vref[i] - Vm))/TA[i];
519 
520     idx = idx + 9;
521   }
522 
523   ierr = VecGetArray(user->V0,&v0);CHKERRQ(ierr);
524   for (i=0; i < nload; i++) {
525     Vr  = xnet[2*lbus[i]]; /* Real part of load bus voltage */
526     Vi  = xnet[2*lbus[i]+1]; /* Imaginary part of the load bus voltage */
527     Vm  = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm;
528     Vm0 = PetscSqrtScalar(v0[2*lbus[i]]*v0[2*lbus[i]] + v0[2*lbus[i]+1]*v0[2*lbus[i]+1]);
529     PD  = QD = 0.0;
530     for (k=0; k < ld_nsegsp[i]; k++) PD += ld_alphap[k]*PD0[i]*PetscPowScalar((Vm/Vm0),ld_betap[k]);
531     for (k=0; k < ld_nsegsq[i]; k++) QD += ld_alphaq[k]*QD0[i]*PetscPowScalar((Vm/Vm0),ld_betaq[k]);
532 
533     /* Load currents */
534     IDr = (PD*Vr + QD*Vi)/Vm2;
535     IDi = (-QD*Vr + PD*Vi)/Vm2;
536 
537     fnet[2*lbus[i]]   += IDi;
538     fnet[2*lbus[i]+1] += IDr;
539   }
540   ierr = VecRestoreArray(user->V0,&v0);CHKERRQ(ierr);
541 
542   ierr = VecRestoreArray(Xgen,&xgen);CHKERRQ(ierr);
543   ierr = VecRestoreArray(Xnet,&xnet);CHKERRQ(ierr);
544   ierr = VecRestoreArray(Fgen,&fgen);CHKERRQ(ierr);
545   ierr = VecRestoreArray(Fnet,&fnet);CHKERRQ(ierr);
546 
547   ierr = DMCompositeGather(user->dmpgrid,INSERT_VALUES,F,Fgen,Fnet);CHKERRQ(ierr);
548   ierr = DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);CHKERRQ(ierr);
549   ierr = DMCompositeRestoreLocalVectors(user->dmpgrid,&Fgen,&Fnet);CHKERRQ(ierr);
550   PetscFunctionReturn(0);
551 }
552 
553 /*   f(x,y)
554      g(x,y)
555  */
556 PetscErrorCode RHSFunction(TS ts,PetscReal t, Vec X, Vec F, void *ctx)
557 {
558   PetscErrorCode ierr;
559   Userctx        *user=(Userctx*)ctx;
560 
561   PetscFunctionBegin;
562   user->t = t;
563   ierr = ResidualFunction(X,F,user);CHKERRQ(ierr);
564   PetscFunctionReturn(0);
565 }
566 
567 /* f(x,y) - \dot{x}
568      g(x,y) = 0
569  */
570 PetscErrorCode IFunction(TS ts,PetscReal t, Vec X, Vec Xdot, Vec F, void *ctx)
571 {
572   PetscErrorCode ierr;
573   PetscScalar    *f,*xdot;
574   PetscInt       i;
575 
576   PetscFunctionBegin;
577 
578   ierr = RHSFunction(ts,t,X,F,ctx);CHKERRQ(ierr);
579   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
580   ierr = VecGetArray(Xdot,&xdot);CHKERRQ(ierr);
581   for (i=0;i < ngen;i++) {
582     f[9*i]   -= xdot[9*i];
583     f[9*i+1] -= xdot[9*i+1];
584     f[9*i+2] -= xdot[9*i+2];
585     f[9*i+3] -= xdot[9*i+3];
586     f[9*i+6] -= xdot[9*i+6];
587     f[9*i+7] -= xdot[9*i+7];
588     f[9*i+8] -= xdot[9*i+8];
589   }
590   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
591   ierr = VecRestoreArray(Xdot,&xdot);CHKERRQ(ierr);
592   PetscFunctionReturn(0);
593 }
594 
595 /* This function is used for solving the algebraic system only during fault on and
596    off times. It computes the entire F and then zeros out the part corresponding to
597    differential equations
598  F = [0;g(y)];
599 */
600 PetscErrorCode AlgFunction(SNES snes, Vec X, Vec F, void *ctx)
601 {
602   PetscErrorCode ierr;
603   Userctx        *user=(Userctx*)ctx;
604   PetscScalar    *f;
605   PetscInt       i;
606 
607   PetscFunctionBegin;
608   ierr = ResidualFunction(X,F,user);CHKERRQ(ierr);
609   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
610   for (i=0; i < ngen; i++) {
611     f[9*i]   = 0;
612     f[9*i+1] = 0;
613     f[9*i+2] = 0;
614     f[9*i+3] = 0;
615     f[9*i+6] = 0;
616     f[9*i+7] = 0;
617     f[9*i+8] = 0;
618   }
619   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
620   PetscFunctionReturn(0);
621 }
622 
623 PetscErrorCode PostStage(TS ts, PetscReal t, PetscInt i, Vec *X)
624 {
625   PetscErrorCode ierr;
626   Userctx        *user;
627 
628   PetscFunctionBegin;
629   ierr = TSGetApplicationContext(ts,&user);CHKERRQ(ierr);
630   ierr = SNESSolve(user->snes_alg,NULL,X[i]);CHKERRQ(ierr);
631   PetscFunctionReturn(0);
632 }
633 
634 PetscErrorCode PostEvaluate(TS ts)
635 {
636   PetscErrorCode ierr;
637   Userctx        *user;
638   Vec            X;
639 
640   PetscFunctionBegin;
641   ierr = TSGetApplicationContext(ts,&user);CHKERRQ(ierr);
642   ierr = TSGetSolution(ts,&X);CHKERRQ(ierr);
643   ierr = SNESSolve(user->snes_alg,NULL,X);CHKERRQ(ierr);
644   PetscFunctionReturn(0);
645 }
646 
647 
648 PetscErrorCode PreallocateJacobian(Mat J, Userctx *user)
649 {
650   PetscErrorCode ierr;
651   PetscInt       *d_nnz;
652   PetscInt       i,idx=0,start=0;
653   PetscInt       ncols;
654 
655   PetscFunctionBegin;
656   ierr = PetscMalloc1(user->neqs_pgrid,&d_nnz);CHKERRQ(ierr);
657   for (i=0; i<user->neqs_pgrid; i++) d_nnz[i] = 0;
658   /* Generator subsystem */
659   for (i=0; i < ngen; i++) {
660 
661     d_nnz[idx]   += 3;
662     d_nnz[idx+1] += 2;
663     d_nnz[idx+2] += 2;
664     d_nnz[idx+3] += 5;
665     d_nnz[idx+4] += 6;
666     d_nnz[idx+5] += 6;
667 
668     d_nnz[user->neqs_gen+2*gbus[i]]   += 3;
669     d_nnz[user->neqs_gen+2*gbus[i]+1] += 3;
670 
671     d_nnz[idx+6] += 2;
672     d_nnz[idx+7] += 2;
673     d_nnz[idx+8] += 5;
674 
675     idx = idx + 9;
676   }
677 
678   start = user->neqs_gen;
679 
680   for (i=0; i < nbus; i++) {
681     ierr = MatGetRow(user->Ybus,2*i,&ncols,NULL,NULL);CHKERRQ(ierr);
682     d_nnz[start+2*i]   += ncols;
683     d_nnz[start+2*i+1] += ncols;
684     ierr = MatRestoreRow(user->Ybus,2*i,&ncols,NULL,NULL);CHKERRQ(ierr);
685   }
686 
687   ierr = MatSeqAIJSetPreallocation(J,0,d_nnz);CHKERRQ(ierr);
688 
689   ierr = PetscFree(d_nnz);CHKERRQ(ierr);
690   PetscFunctionReturn(0);
691 }
692 
693 /*
694    J = [df_dx, df_dy
695         dg_dx, dg_dy]
696 */
697 PetscErrorCode ResidualJacobian(Vec X,Mat J,Mat B,void *ctx)
698 {
699   PetscErrorCode ierr;
700   Userctx        *user=(Userctx*)ctx;
701   Vec            Xgen,Xnet;
702   PetscScalar    *xgen,*xnet;
703   PetscInt       i,idx=0;
704   PetscScalar    Vr,Vi,Vm,Vm2;
705   PetscScalar    Eqp,Edp,delta; /* Generator variables */
706   PetscScalar    Efd;
707   PetscScalar    Id,Iq;  /* Generator dq axis currents */
708   PetscScalar    Vd,Vq;
709   PetscScalar    val[10];
710   PetscInt       row[2],col[10];
711   PetscInt       net_start=user->neqs_gen;
712   PetscScalar    Zdq_inv[4],det;
713   PetscScalar    dVd_dVr,dVd_dVi,dVq_dVr,dVq_dVi,dVd_ddelta,dVq_ddelta;
714   PetscScalar    dIGr_ddelta,dIGi_ddelta,dIGr_dId,dIGr_dIq,dIGi_dId,dIGi_dIq;
715   PetscScalar    dSE_dEfd;
716   PetscScalar    dVm_dVd,dVm_dVq,dVm_dVr,dVm_dVi;
717   PetscInt          ncols;
718   const PetscInt    *cols;
719   const PetscScalar *yvals;
720   PetscInt          k;
721   PetscScalar PD,QD,Vm0,*v0,Vm4;
722   PetscScalar dPD_dVr,dPD_dVi,dQD_dVr,dQD_dVi;
723   PetscScalar dIDr_dVr,dIDr_dVi,dIDi_dVr,dIDi_dVi;
724 
725 
726   PetscFunctionBegin;
727   ierr  = MatZeroEntries(B);CHKERRQ(ierr);
728   ierr  = DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);CHKERRQ(ierr);
729   ierr  = DMCompositeScatter(user->dmpgrid,X,Xgen,Xnet);CHKERRQ(ierr);
730 
731   ierr = VecGetArray(Xgen,&xgen);CHKERRQ(ierr);
732   ierr = VecGetArray(Xnet,&xnet);CHKERRQ(ierr);
733 
734   /* Generator subsystem */
735   for (i=0; i < ngen; i++) {
736     Eqp   = xgen[idx];
737     Edp   = xgen[idx+1];
738     delta = xgen[idx+2];
739     Id    = xgen[idx+4];
740     Iq    = xgen[idx+5];
741     Efd   = xgen[idx+6];
742 
743     /*    fgen[idx]   = (-Eqp - (Xd[i] - Xdp[i])*Id + Efd)/Td0p[i]; */
744     row[0] = idx;
745     col[0] = idx;           col[1] = idx+4;          col[2] = idx+6;
746     val[0] = -1/ Td0p[i]; val[1] = -(Xd[i] - Xdp[i])/ Td0p[i]; val[2] = 1/Td0p[i];
747 
748     ierr = MatSetValues(J,1,row,3,col,val,INSERT_VALUES);CHKERRQ(ierr);
749 
750     /*    fgen[idx+1] = (-Edp + (Xq[i] - Xqp[i])*Iq)/Tq0p[i]; */
751     row[0] = idx + 1;
752     col[0] = idx + 1;       col[1] = idx+5;
753     val[0] = -1/Tq0p[i]; val[1] = (Xq[i] - Xqp[i])/Tq0p[i];
754     ierr   = MatSetValues(J,1,row,2,col,val,INSERT_VALUES);CHKERRQ(ierr);
755 
756     /*    fgen[idx+2] = w - w_s; */
757     row[0] = idx + 2;
758     col[0] = idx + 2; col[1] = idx + 3;
759     val[0] = 0;       val[1] = 1;
760     ierr   = MatSetValues(J,1,row,2,col,val,INSERT_VALUES);CHKERRQ(ierr);
761 
762     /*    fgen[idx+3] = (TM[i] - Edp*Id - Eqp*Iq - (Xqp[i] - Xdp[i])*Id*Iq - D[i]*(w - w_s))/M[i]; */
763     row[0] = idx + 3;
764     col[0] = idx; col[1] = idx + 1; col[2] = idx + 3;       col[3] = idx + 4;                  col[4] = idx + 5;
765     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];
766     ierr   = MatSetValues(J,1,row,5,col,val,INSERT_VALUES);CHKERRQ(ierr);
767 
768     Vr   = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */
769     Vi   = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */
770     ierr = ri2dq(Vr,Vi,delta,&Vd,&Vq);CHKERRQ(ierr);
771 
772     det = Rs[i]*Rs[i] + Xdp[i]*Xqp[i];
773 
774     Zdq_inv[0] = Rs[i]/det;
775     Zdq_inv[1] = Xqp[i]/det;
776     Zdq_inv[2] = -Xdp[i]/det;
777     Zdq_inv[3] = Rs[i]/det;
778 
779     dVd_dVr    = PetscSinScalar(delta); dVd_dVi = -PetscCosScalar(delta);
780     dVq_dVr    = PetscCosScalar(delta); dVq_dVi = PetscSinScalar(delta);
781     dVd_ddelta = Vr*PetscCosScalar(delta) + Vi*PetscSinScalar(delta);
782     dVq_ddelta = -Vr*PetscSinScalar(delta) + Vi*PetscCosScalar(delta);
783 
784     /*    fgen[idx+4] = Zdq_inv[0]*(-Edp + Vd) + Zdq_inv[1]*(-Eqp + Vq) + Id; */
785     row[0] = idx+4;
786     col[0] = idx;         col[1] = idx+1;        col[2] = idx + 2;
787     val[0] = -Zdq_inv[1]; val[1] = -Zdq_inv[0];  val[2] = Zdq_inv[0]*dVd_ddelta + Zdq_inv[1]*dVq_ddelta;
788     col[3] = idx + 4; col[4] = net_start+2*gbus[i];                     col[5] = net_start + 2*gbus[i]+1;
789     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;
790     ierr   = MatSetValues(J,1,row,6,col,val,INSERT_VALUES);CHKERRQ(ierr);
791 
792     /*  fgen[idx+5] = Zdq_inv[2]*(-Edp + Vd) + Zdq_inv[3]*(-Eqp + Vq) + Iq; */
793     row[0] = idx+5;
794     col[0] = idx;         col[1] = idx+1;        col[2] = idx + 2;
795     val[0] = -Zdq_inv[3]; val[1] = -Zdq_inv[2];  val[2] = Zdq_inv[2]*dVd_ddelta + Zdq_inv[3]*dVq_ddelta;
796     col[3] = idx + 5; col[4] = net_start+2*gbus[i];                     col[5] = net_start + 2*gbus[i]+1;
797     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;
798     ierr   = MatSetValues(J,1,row,6,col,val,INSERT_VALUES);CHKERRQ(ierr);
799 
800     dIGr_ddelta = Id*PetscCosScalar(delta) - Iq*PetscSinScalar(delta);
801     dIGi_ddelta = Id*PetscSinScalar(delta) + Iq*PetscCosScalar(delta);
802     dIGr_dId    = PetscSinScalar(delta);  dIGr_dIq = PetscCosScalar(delta);
803     dIGi_dId    = -PetscCosScalar(delta); dIGi_dIq = PetscSinScalar(delta);
804 
805     /* fnet[2*gbus[i]]   -= IGi; */
806     row[0] = net_start + 2*gbus[i];
807     col[0] = idx+2;        col[1] = idx + 4;   col[2] = idx + 5;
808     val[0] = -dIGi_ddelta; val[1] = -dIGi_dId; val[2] = -dIGi_dIq;
809     ierr = MatSetValues(J,1,row,3,col,val,INSERT_VALUES);CHKERRQ(ierr);
810 
811     /* fnet[2*gbus[i]+1]   -= IGr; */
812     row[0] = net_start + 2*gbus[i]+1;
813     col[0] = idx+2;        col[1] = idx + 4;   col[2] = idx + 5;
814     val[0] = -dIGr_ddelta; val[1] = -dIGr_dId; val[2] = -dIGr_dIq;
815     ierr   = MatSetValues(J,1,row,3,col,val,INSERT_VALUES);CHKERRQ(ierr);
816 
817     Vm = PetscSqrtScalar(Vd*Vd + Vq*Vq);
818 
819     /*    fgen[idx+6] = (-KE[i]*Efd - SE + VR)/TE[i]; */
820     /*    SE  = k1[i]*PetscExpScalar(k2[i]*Efd); */
821 
822     dSE_dEfd = k1[i]*k2[i]*PetscExpScalar(k2[i]*Efd);
823 
824     row[0] = idx + 6;
825     col[0] = idx + 6;                     col[1] = idx + 8;
826     val[0] = (-KE[i] - dSE_dEfd)/TE[i];  val[1] = 1/TE[i];
827     ierr   = MatSetValues(J,1,row,2,col,val,INSERT_VALUES);CHKERRQ(ierr);
828 
829     /* Exciter differential equations */
830 
831     /*    fgen[idx+7] = (-RF + KF[i]*Efd/TF[i])/TF[i]; */
832     row[0] = idx + 7;
833     col[0] = idx + 6;       col[1] = idx + 7;
834     val[0] = (KF[i]/TF[i])/TF[i];  val[1] = -1/TF[i];
835     ierr   = MatSetValues(J,1,row,2,col,val,INSERT_VALUES);CHKERRQ(ierr);
836 
837     /*    fgen[idx+8] = (-VR + KA[i]*RF - KA[i]*KF[i]*Efd/TF[i] + KA[i]*(Vref[i] - Vm))/TA[i]; */
838     /* Vm = (Vd^2 + Vq^2)^0.5; */
839 
840     row[0] = idx + 8;
841     if(VRatmax[i]) {
842       col[0] = idx + 8; val[0] = 1.0;
843       ierr = MatSetValues(J,1,row,1,col,val,INSERT_VALUES);CHKERRQ(ierr);
844     } else if(VRatmin[i]) {
845       col[0] = idx + 8; val[0] = -1.0;
846       ierr = MatSetValues(J,1,row,1,col,val,INSERT_VALUES);CHKERRQ(ierr);
847     } else {
848       dVm_dVd    = Vd/Vm; dVm_dVq = Vq/Vm;
849       dVm_dVr    = dVm_dVd*dVd_dVr + dVm_dVq*dVq_dVr;
850       dVm_dVi    = dVm_dVd*dVd_dVi + dVm_dVq*dVq_dVi;
851       row[0]     = idx + 8;
852       col[0]     = idx + 6;           col[1] = idx + 7; col[2] = idx + 8;
853       val[0]     = -(KA[i]*KF[i]/TF[i])/TA[i]; val[1] = KA[i]/TA[i];  val[2] = -1/TA[i];
854       col[3]     = net_start + 2*gbus[i]; col[4] = net_start + 2*gbus[i]+1;
855       val[3]     = -KA[i]*dVm_dVr/TA[i];         val[4] = -KA[i]*dVm_dVi/TA[i];
856       ierr       = MatSetValues(J,1,row,5,col,val,INSERT_VALUES);CHKERRQ(ierr);
857     }
858 
859     idx        = idx + 9;
860   }
861 
862   for (i=0; i<nbus; i++) {
863     ierr   = MatGetRow(user->Ybus,2*i,&ncols,&cols,&yvals);CHKERRQ(ierr);
864     row[0] = net_start + 2*i;
865     for (k=0; k<ncols; k++) {
866       col[k] = net_start + cols[k];
867       val[k] = yvals[k];
868     }
869     ierr = MatSetValues(J,1,row,ncols,col,val,INSERT_VALUES);CHKERRQ(ierr);
870     ierr = MatRestoreRow(user->Ybus,2*i,&ncols,&cols,&yvals);CHKERRQ(ierr);
871 
872     ierr   = MatGetRow(user->Ybus,2*i+1,&ncols,&cols,&yvals);CHKERRQ(ierr);
873     row[0] = net_start + 2*i+1;
874     for (k=0; k<ncols; k++) {
875       col[k] = net_start + cols[k];
876       val[k] = yvals[k];
877     }
878     ierr = MatSetValues(J,1,row,ncols,col,val,INSERT_VALUES);CHKERRQ(ierr);
879     ierr = MatRestoreRow(user->Ybus,2*i+1,&ncols,&cols,&yvals);CHKERRQ(ierr);
880   }
881 
882   ierr = MatAssemblyBegin(J,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr);
883   ierr = MatAssemblyEnd(J,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr);
884 
885   ierr = VecGetArray(user->V0,&v0);CHKERRQ(ierr);
886   for (i=0; i < nload; i++) {
887     Vr      = xnet[2*lbus[i]]; /* Real part of load bus voltage */
888     Vi      = xnet[2*lbus[i]+1]; /* Imaginary part of the load bus voltage */
889     Vm      = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm; Vm4 = Vm2*Vm2;
890     Vm0     = PetscSqrtScalar(v0[2*lbus[i]]*v0[2*lbus[i]] + v0[2*lbus[i]+1]*v0[2*lbus[i]+1]);
891     PD      = QD = 0.0;
892     dPD_dVr = dPD_dVi = dQD_dVr = dQD_dVi = 0.0;
893     for (k=0; k < ld_nsegsp[i]; k++) {
894       PD      += ld_alphap[k]*PD0[i]*PetscPowScalar((Vm/Vm0),ld_betap[k]);
895       dPD_dVr += ld_alphap[k]*ld_betap[k]*PD0[i]*PetscPowScalar((1/Vm0),ld_betap[k])*Vr*PetscPowScalar(Vm,(ld_betap[k]-2));
896       dPD_dVi += ld_alphap[k]*ld_betap[k]*PD0[i]*PetscPowScalar((1/Vm0),ld_betap[k])*Vi*PetscPowScalar(Vm,(ld_betap[k]-2));
897     }
898     for (k=0; k < ld_nsegsq[i]; k++) {
899       QD      += ld_alphaq[k]*QD0[i]*PetscPowScalar((Vm/Vm0),ld_betaq[k]);
900       dQD_dVr += ld_alphaq[k]*ld_betaq[k]*QD0[i]*PetscPowScalar((1/Vm0),ld_betaq[k])*Vr*PetscPowScalar(Vm,(ld_betaq[k]-2));
901       dQD_dVi += ld_alphaq[k]*ld_betaq[k]*QD0[i]*PetscPowScalar((1/Vm0),ld_betaq[k])*Vi*PetscPowScalar(Vm,(ld_betaq[k]-2));
902     }
903 
904     /*    IDr = (PD*Vr + QD*Vi)/Vm2; */
905     /*    IDi = (-QD*Vr + PD*Vi)/Vm2; */
906 
907     dIDr_dVr = (dPD_dVr*Vr + dQD_dVr*Vi + PD)/Vm2 - ((PD*Vr + QD*Vi)*2*Vr)/Vm4;
908     dIDr_dVi = (dPD_dVi*Vr + dQD_dVi*Vi + QD)/Vm2 - ((PD*Vr + QD*Vi)*2*Vi)/Vm4;
909 
910     dIDi_dVr = (-dQD_dVr*Vr + dPD_dVr*Vi - QD)/Vm2 - ((-QD*Vr + PD*Vi)*2*Vr)/Vm4;
911     dIDi_dVi = (-dQD_dVi*Vr + dPD_dVi*Vi + PD)/Vm2 - ((-QD*Vr + PD*Vi)*2*Vi)/Vm4;
912 
913 
914     /*    fnet[2*lbus[i]]   += IDi; */
915     row[0] = net_start + 2*lbus[i];
916     col[0] = net_start + 2*lbus[i];  col[1] = net_start + 2*lbus[i]+1;
917     val[0] = dIDi_dVr;               val[1] = dIDi_dVi;
918     ierr   = MatSetValues(J,1,row,2,col,val,ADD_VALUES);CHKERRQ(ierr);
919     /*    fnet[2*lbus[i]+1] += IDr; */
920     row[0] = net_start + 2*lbus[i]+1;
921     col[0] = net_start + 2*lbus[i];  col[1] = net_start + 2*lbus[i]+1;
922     val[0] = dIDr_dVr;               val[1] = dIDr_dVi;
923     ierr   = MatSetValues(J,1,row,2,col,val,ADD_VALUES);CHKERRQ(ierr);
924   }
925   ierr = VecRestoreArray(user->V0,&v0);CHKERRQ(ierr);
926 
927   ierr = VecRestoreArray(Xgen,&xgen);CHKERRQ(ierr);
928   ierr = VecRestoreArray(Xnet,&xnet);CHKERRQ(ierr);
929 
930   ierr = DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);CHKERRQ(ierr);
931 
932   ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
933   ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
934   PetscFunctionReturn(0);
935 }
936 
937 /*
938    J = [I, 0
939         dg_dx, dg_dy]
940 */
941 PetscErrorCode AlgJacobian(SNES snes,Vec X,Mat A,Mat B,void *ctx)
942 {
943   PetscErrorCode ierr;
944   Userctx        *user=(Userctx*)ctx;
945 
946   PetscFunctionBegin;
947   ierr = ResidualJacobian(X,A,B,ctx);CHKERRQ(ierr);
948   ierr = MatSetOption(A,MAT_KEEP_NONZERO_PATTERN,PETSC_TRUE);CHKERRQ(ierr);
949   ierr = MatZeroRowsIS(A,user->is_diff,1.0,NULL,NULL);CHKERRQ(ierr);
950   PetscFunctionReturn(0);
951 }
952 
953 /*
954    J = [-df_dx, -df_dy
955         dg_dx, dg_dy]
956 */
957 
958 PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec X,Mat A,Mat B,void *ctx)
959 {
960   PetscErrorCode ierr;
961   Userctx        *user=(Userctx*)ctx;
962 
963   PetscFunctionBegin;
964   user->t = t;
965 
966   ierr = ResidualJacobian(X,A,B,user);CHKERRQ(ierr);
967 
968   PetscFunctionReturn(0);
969 }
970 
971 /*
972    J = [df_dx-aI, df_dy
973         dg_dx, dg_dy]
974 */
975 
976 PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,Userctx *user)
977 {
978   PetscErrorCode ierr;
979   PetscScalar    atmp = (PetscScalar) a;
980   PetscInt       i,row;
981 
982   PetscFunctionBegin;
983   user->t = t;
984   atmp *= -1;
985 
986   ierr = RHSJacobian(ts,t,X,A,B,user);CHKERRQ(ierr);
987   for (i=0;i < ngen;i++) {
988     row = 9*i;
989     ierr = MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);CHKERRQ(ierr);
990     row  = 9*i+1;
991     ierr = MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);CHKERRQ(ierr);
992     row  = 9*i+2;
993     ierr = MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);CHKERRQ(ierr);
994     row  = 9*i+3;
995     ierr = MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);CHKERRQ(ierr);
996     row  = 9*i+6;
997     ierr = MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);CHKERRQ(ierr);
998     row  = 9*i+7;
999     ierr = MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);CHKERRQ(ierr);
1000     row  = 9*i+8;
1001     ierr = MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);CHKERRQ(ierr);
1002   }
1003   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1004   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1005   PetscFunctionReturn(0);
1006 }
1007 
1008 int main(int argc,char **argv)
1009 {
1010   TS               ts;
1011   SNES              snes_alg;
1012   PetscErrorCode    ierr;
1013   PetscMPIInt       size;
1014   Userctx           user;
1015   PetscViewer       Xview,Ybusview,viewer;
1016   Vec               X,F_alg;
1017   Mat               J,A;
1018   PetscInt          i,idx,*idx2;
1019   Vec               Xdot;
1020   PetscScalar       *x,*mat,*amat;
1021   const PetscScalar *rmat;
1022   Vec               vatol;
1023   PetscInt          *direction;
1024   PetscBool         *terminate;
1025   const PetscInt    *idx3;
1026   PetscScalar       *vatoli;
1027   PetscInt          k;
1028 
1029 
1030   ierr = PetscInitialize(&argc,&argv,"petscoptions",help);if (ierr) return ierr;
1031   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
1032   if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs");
1033 
1034   user.neqs_gen   = 9*ngen; /* # eqs. for generator subsystem */
1035   user.neqs_net   = 2*nbus; /* # eqs. for network subsystem   */
1036   user.neqs_pgrid = user.neqs_gen + user.neqs_net;
1037 
1038   /* Create indices for differential and algebraic equations */
1039 
1040   ierr = PetscMalloc1(7*ngen,&idx2);CHKERRQ(ierr);
1041   for (i=0; i<ngen; i++) {
1042     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;
1043     idx2[7*i+4] = 9*i+6; idx2[7*i+5] = 9*i+7; idx2[7*i+6] = 9*i+8;
1044   }
1045   ierr = ISCreateGeneral(PETSC_COMM_WORLD,7*ngen,idx2,PETSC_COPY_VALUES,&user.is_diff);CHKERRQ(ierr);
1046   ierr = ISComplement(user.is_diff,0,user.neqs_pgrid,&user.is_alg);CHKERRQ(ierr);
1047   ierr = PetscFree(idx2);CHKERRQ(ierr);
1048 
1049   /* Read initial voltage vector and Ybus */
1050   ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,"X.bin",FILE_MODE_READ,&Xview);CHKERRQ(ierr);
1051   ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,"Ybus.bin",FILE_MODE_READ,&Ybusview);CHKERRQ(ierr);
1052 
1053   ierr = VecCreate(PETSC_COMM_WORLD,&user.V0);CHKERRQ(ierr);
1054   ierr = VecSetSizes(user.V0,PETSC_DECIDE,user.neqs_net);CHKERRQ(ierr);
1055   ierr = VecLoad(user.V0,Xview);CHKERRQ(ierr);
1056 
1057   ierr = MatCreate(PETSC_COMM_WORLD,&user.Ybus);CHKERRQ(ierr);
1058   ierr = MatSetSizes(user.Ybus,PETSC_DECIDE,PETSC_DECIDE,user.neqs_net,user.neqs_net);CHKERRQ(ierr);
1059   ierr = MatSetType(user.Ybus,MATBAIJ);CHKERRQ(ierr);
1060   /*  ierr = MatSetBlockSize(user.Ybus,2);CHKERRQ(ierr); */
1061   ierr = MatLoad(user.Ybus,Ybusview);CHKERRQ(ierr);
1062 
1063   /* Set run time options */
1064   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Transient stability fault options","");CHKERRQ(ierr);
1065   {
1066     user.tfaulton  = 1.0;
1067     user.tfaultoff = 1.2;
1068     user.Rfault    = 0.0001;
1069     user.setisdiff = PETSC_FALSE;
1070     user.semiexplicit = PETSC_FALSE;
1071     user.faultbus  = 8;
1072     ierr           = PetscOptionsReal("-tfaulton","","",user.tfaulton,&user.tfaulton,NULL);CHKERRQ(ierr);
1073     ierr           = PetscOptionsReal("-tfaultoff","","",user.tfaultoff,&user.tfaultoff,NULL);CHKERRQ(ierr);
1074     ierr           = PetscOptionsInt("-faultbus","","",user.faultbus,&user.faultbus,NULL);CHKERRQ(ierr);
1075     user.t0        = 0.0;
1076     user.tmax      = 5.0;
1077     ierr           = PetscOptionsReal("-t0","","",user.t0,&user.t0,NULL);CHKERRQ(ierr);
1078     ierr           = PetscOptionsReal("-tmax","","",user.tmax,&user.tmax,NULL);CHKERRQ(ierr);
1079     ierr           = PetscOptionsBool("-setisdiff","","",user.setisdiff,&user.setisdiff,NULL);CHKERRQ(ierr);
1080     ierr           = PetscOptionsBool("-dae_semiexplicit","","",user.semiexplicit,&user.semiexplicit,NULL);CHKERRQ(ierr);
1081   }
1082   ierr = PetscOptionsEnd();CHKERRQ(ierr);
1083 
1084   ierr = PetscViewerDestroy(&Xview);CHKERRQ(ierr);
1085   ierr = PetscViewerDestroy(&Ybusview);CHKERRQ(ierr);
1086 
1087   /* Create DMs for generator and network subsystems */
1088   ierr = DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,user.neqs_gen,1,1,NULL,&user.dmgen);CHKERRQ(ierr);
1089   ierr = DMSetOptionsPrefix(user.dmgen,"dmgen_");CHKERRQ(ierr);
1090   ierr = DMSetFromOptions(user.dmgen);CHKERRQ(ierr);
1091   ierr = DMSetUp(user.dmgen);CHKERRQ(ierr);
1092   ierr = DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,user.neqs_net,1,1,NULL,&user.dmnet);CHKERRQ(ierr);
1093   ierr = DMSetOptionsPrefix(user.dmnet,"dmnet_");CHKERRQ(ierr);
1094   ierr = DMSetFromOptions(user.dmnet);CHKERRQ(ierr);
1095   ierr = DMSetUp(user.dmnet);CHKERRQ(ierr);
1096   /* Create a composite DM packer and add the two DMs */
1097   ierr = DMCompositeCreate(PETSC_COMM_WORLD,&user.dmpgrid);CHKERRQ(ierr);
1098   ierr = DMSetOptionsPrefix(user.dmpgrid,"pgrid_");CHKERRQ(ierr);
1099   ierr = DMCompositeAddDM(user.dmpgrid,user.dmgen);CHKERRQ(ierr);
1100   ierr = DMCompositeAddDM(user.dmpgrid,user.dmnet);CHKERRQ(ierr);
1101 
1102   ierr = DMCreateGlobalVector(user.dmpgrid,&X);CHKERRQ(ierr);
1103 
1104   ierr = MatCreate(PETSC_COMM_WORLD,&J);CHKERRQ(ierr);
1105   ierr = MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,user.neqs_pgrid,user.neqs_pgrid);CHKERRQ(ierr);
1106   ierr = MatSetFromOptions(J);CHKERRQ(ierr);
1107   ierr = PreallocateJacobian(J,&user);CHKERRQ(ierr);
1108 
1109   /* Create matrix to save solutions at each time step */
1110   user.stepnum = 0;
1111 
1112   ierr = MatCreateSeqDense(PETSC_COMM_SELF,user.neqs_pgrid+1,1002,NULL,&user.Sol);CHKERRQ(ierr);
1113   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1114      Create timestepping solver context
1115      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1116   ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
1117   ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr);
1118   if(user.semiexplicit) {
1119     ierr = TSSetType(ts,TSRK);CHKERRQ(ierr);
1120     ierr = TSSetRHSFunction(ts,NULL,RHSFunction,&user);CHKERRQ(ierr);
1121     ierr = TSSetRHSJacobian(ts,J,J,RHSJacobian,&user);CHKERRQ(ierr);
1122   } else {
1123     ierr = TSSetType(ts,TSCN);CHKERRQ(ierr);
1124     ierr = TSSetEquationType(ts,TS_EQ_DAE_IMPLICIT_INDEX1);CHKERRQ(ierr);
1125     ierr = TSARKIMEXSetFullyImplicit(ts,PETSC_TRUE);CHKERRQ(ierr);
1126     ierr = TSSetIFunction(ts,NULL,(TSIFunction) IFunction,&user);CHKERRQ(ierr);
1127     ierr = TSSetIJacobian(ts,J,J,(TSIJacobian)IJacobian,&user);CHKERRQ(ierr);
1128   }
1129   ierr = TSSetApplicationContext(ts,&user);CHKERRQ(ierr);
1130 
1131   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1132      Set initial conditions
1133    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1134   ierr = SetInitialGuess(X,&user);CHKERRQ(ierr);
1135   /* Just to set up the Jacobian structure */
1136 
1137   ierr = VecDuplicate(X,&Xdot);CHKERRQ(ierr);
1138   ierr = IJacobian(ts,0.0,X,Xdot,0.0,J,J,&user);CHKERRQ(ierr);
1139   ierr = VecDestroy(&Xdot);CHKERRQ(ierr);
1140 
1141   /* Save initial solution */
1142 
1143   idx=user.stepnum*(user.neqs_pgrid+1);
1144   ierr = MatDenseGetArray(user.Sol,&mat);CHKERRQ(ierr);
1145   ierr = VecGetArray(X,&x);CHKERRQ(ierr);
1146 
1147   mat[idx] = 0.0;
1148 
1149   ierr = PetscArraycpy(mat+idx+1,x,user.neqs_pgrid);CHKERRQ(ierr);
1150   ierr = MatDenseRestoreArray(user.Sol,&mat);CHKERRQ(ierr);
1151   ierr = VecRestoreArray(X,&x);CHKERRQ(ierr);
1152   user.stepnum++;
1153 
1154   ierr = TSSetMaxTime(ts,user.tmax);CHKERRQ(ierr);
1155   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr);
1156   ierr = TSSetTimeStep(ts,0.01);CHKERRQ(ierr);
1157   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
1158   ierr = TSSetPostStep(ts,SaveSolution);CHKERRQ(ierr);
1159   ierr = TSSetSolution(ts,X);CHKERRQ(ierr);
1160 
1161   ierr = PetscMalloc1((2*ngen+2),&direction);CHKERRQ(ierr);
1162   ierr = PetscMalloc1((2*ngen+2),&terminate);CHKERRQ(ierr);
1163   direction[0] = direction[1] = 1;
1164   terminate[0] = terminate[1] = PETSC_FALSE;
1165   for (i=0; i < ngen;i++) {
1166     direction[2+2*i] = -1; direction[2+2*i+1] = 1;
1167     terminate[2+2*i] = terminate[2+2*i+1] = PETSC_FALSE;
1168   }
1169 
1170   ierr = TSSetEventHandler(ts,2*ngen+2,direction,terminate,EventFunction,PostEventFunction,(void*)&user);CHKERRQ(ierr);
1171 
1172   if(user.semiexplicit) {
1173     /* Use a semi-explicit approach with the time-stepping done by an explicit method and the
1174        algrebraic part solved via PostStage and PostEvaluate callbacks
1175     */
1176     ierr = TSSetType(ts,TSRK);CHKERRQ(ierr);
1177     ierr = TSSetPostStage(ts,PostStage);CHKERRQ(ierr);
1178     ierr = TSSetPostEvaluate(ts,PostEvaluate);CHKERRQ(ierr);
1179   }
1180 
1181 
1182   if(user.setisdiff) {
1183     /* Create vector of absolute tolerances and set the algebraic part to infinity */
1184     ierr = VecDuplicate(X,&vatol);CHKERRQ(ierr);
1185     ierr = VecSet(vatol,100000.0);CHKERRQ(ierr);
1186     ierr = VecGetArray(vatol,&vatoli);CHKERRQ(ierr);
1187     ierr = ISGetIndices(user.is_diff,&idx3);CHKERRQ(ierr);
1188     for(k=0; k < 7*ngen; k++) vatoli[idx3[k]] = 1e-2;
1189     ierr = VecRestoreArray(vatol,&vatoli);CHKERRQ(ierr);
1190   }
1191 
1192   /* Create the nonlinear solver for solving the algebraic system */
1193   /* Note that although the algebraic system needs to be solved only for
1194      Idq and V, we reuse the entire system including xgen. The xgen
1195      variables are held constant by setting their residuals to 0 and
1196      putting a 1 on the Jacobian diagonal for xgen rows
1197   */
1198 
1199   ierr = VecDuplicate(X,&F_alg);CHKERRQ(ierr);
1200   ierr = SNESCreate(PETSC_COMM_WORLD,&snes_alg);CHKERRQ(ierr);
1201   ierr = SNESSetFunction(snes_alg,F_alg,AlgFunction,&user);CHKERRQ(ierr);
1202   ierr = SNESSetJacobian(snes_alg,J,J,AlgJacobian,&user);CHKERRQ(ierr);
1203   ierr = SNESSetFromOptions(snes_alg);CHKERRQ(ierr);
1204 
1205   user.snes_alg=snes_alg;
1206 
1207   /* Solve */
1208   ierr = TSSolve(ts,X);CHKERRQ(ierr);
1209 
1210   ierr = MatAssemblyBegin(user.Sol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1211   ierr = MatAssemblyEnd(user.Sol,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1212 
1213   ierr = MatCreateSeqDense(PETSC_COMM_SELF,user.neqs_pgrid+1,user.stepnum,NULL,&A);CHKERRQ(ierr);
1214   ierr = MatDenseGetArrayRead(user.Sol,&rmat);CHKERRQ(ierr);
1215   ierr = MatDenseGetArray(A,&amat);CHKERRQ(ierr);
1216   ierr = PetscArraycpy(amat,rmat,user.stepnum*(user.neqs_pgrid+1));CHKERRQ(ierr);
1217   ierr = MatDenseRestoreArray(A,&amat);CHKERRQ(ierr);
1218   ierr = MatDenseRestoreArrayRead(user.Sol,&rmat);CHKERRQ(ierr);
1219   ierr = PetscViewerBinaryOpen(PETSC_COMM_SELF,"out.bin",FILE_MODE_WRITE,&viewer);CHKERRQ(ierr);
1220   ierr = MatView(A,viewer);CHKERRQ(ierr);
1221   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1222   ierr = MatDestroy(&A);CHKERRQ(ierr);
1223 
1224   ierr = PetscFree(direction);CHKERRQ(ierr);
1225   ierr = PetscFree(terminate);CHKERRQ(ierr);
1226   ierr = SNESDestroy(&snes_alg);CHKERRQ(ierr);
1227   ierr = VecDestroy(&F_alg);CHKERRQ(ierr);
1228   ierr = MatDestroy(&J);CHKERRQ(ierr);
1229   ierr = MatDestroy(&user.Ybus);CHKERRQ(ierr);
1230   ierr = MatDestroy(&user.Sol);CHKERRQ(ierr);
1231   ierr = VecDestroy(&X);CHKERRQ(ierr);
1232   ierr = VecDestroy(&user.V0);CHKERRQ(ierr);
1233   ierr = DMDestroy(&user.dmgen);CHKERRQ(ierr);
1234   ierr = DMDestroy(&user.dmnet);CHKERRQ(ierr);
1235   ierr = DMDestroy(&user.dmpgrid);CHKERRQ(ierr);
1236   ierr = ISDestroy(&user.is_diff);CHKERRQ(ierr);
1237   ierr = ISDestroy(&user.is_alg);CHKERRQ(ierr);
1238   ierr = TSDestroy(&ts);CHKERRQ(ierr);
1239   if(user.setisdiff) {
1240     ierr = VecDestroy(&vatol);CHKERRQ(ierr);
1241   }
1242   ierr = PetscFinalize();
1243   return ierr;
1244 }
1245 
1246 /*TEST
1247 
1248    build:
1249       requires: double !complex !define(PETSC_USE_64BIT_INDICES)
1250 
1251    test:
1252       suffix: implicit
1253       args: -ts_monitor -snes_monitor_short
1254       localrunfiles: petscoptions X.bin Ybus.bin
1255 
1256    test:
1257       suffix: semiexplicit
1258       args: -ts_monitor -snes_monitor_short -dae_semiexplicit -ts_rk_type 2a
1259       localrunfiles: petscoptions X.bin Ybus.bin
1260 
1261    test:
1262       suffix: steprestart
1263       args: -ts_monitor -snes_monitor_short -ts_type arkimex
1264       localrunfiles: petscoptions X.bin Ybus.bin
1265 
1266 TEST*/
1267