xref: /petsc/src/ts/tutorials/power_grid/stability_9bus/ex9bus.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
1c4762a1bSJed Brown static char help[] = "Power grid stability analysis of WECC 9 bus system.\n\
2c4762a1bSJed Brown This example is based on the 9-bus (node) example given in the book Power\n\
3c4762a1bSJed Brown Systems Dynamics and Stability (Chapter 7) by P. Sauer and M. A. Pai.\n\
4c4762a1bSJed Brown The power grid in this example consists of 9 buses (nodes), 3 generators,\n\
5c4762a1bSJed Brown 3 loads, and 9 transmission lines. The network equations are written\n\
6a5b23f4aSJose E. Roman in current balance form using rectangular coordinates.\n\n";
7c4762a1bSJed Brown 
8c4762a1bSJed Brown /*
9c4762a1bSJed Brown    The equations for the stability analysis are described by the DAE
10c4762a1bSJed Brown 
11c4762a1bSJed Brown    \dot{x} = f(x,y,t)
12c4762a1bSJed Brown      0     = g(x,y,t)
13c4762a1bSJed Brown 
14c4762a1bSJed Brown    where the generators are described by differential equations, while the algebraic
15c4762a1bSJed Brown    constraints define the network equations.
16c4762a1bSJed Brown 
17c4762a1bSJed Brown    The generators are modeled with a 4th order differential equation describing the electrical
18c4762a1bSJed Brown    and mechanical dynamics. Each generator also has an exciter system modeled by 3rd order
19c4762a1bSJed Brown    diff. eqns. describing the exciter, voltage regulator, and the feedback stabilizer
20c4762a1bSJed Brown    mechanism.
21c4762a1bSJed Brown 
22c4762a1bSJed Brown    The network equations are described by nodal current balance equations.
23c4762a1bSJed Brown     I(x,y) - Y*V = 0
24c4762a1bSJed Brown 
25c4762a1bSJed Brown    where:
26c4762a1bSJed Brown     I(x,y) is the current injected from generators and loads.
27c4762a1bSJed Brown       Y    is the admittance matrix, and
28c4762a1bSJed Brown       V    is the voltage vector
29c4762a1bSJed Brown */
30c4762a1bSJed Brown 
31c4762a1bSJed Brown /*
32c4762a1bSJed Brown    Include "petscts.h" so that we can use TS solvers.  Note that this
33c4762a1bSJed Brown    file automatically includes:
34c4762a1bSJed Brown      petscsys.h       - base PETSc routines   petscvec.h - vectors
35c4762a1bSJed Brown      petscmat.h - matrices
36c4762a1bSJed Brown      petscis.h     - index sets            petscksp.h - Krylov subspace methods
37c4762a1bSJed Brown      petscviewer.h - viewers               petscpc.h  - preconditioners
38c4762a1bSJed Brown      petscksp.h   - linear solvers
39c4762a1bSJed Brown */
40c4762a1bSJed Brown 
41c4762a1bSJed Brown #include <petscts.h>
42c4762a1bSJed Brown #include <petscdm.h>
43c4762a1bSJed Brown #include <petscdmda.h>
44c4762a1bSJed Brown #include <petscdmcomposite.h>
45c4762a1bSJed Brown 
46c4762a1bSJed Brown #define freq 60
47c4762a1bSJed Brown #define w_s  (2 * PETSC_PI * freq)
48c4762a1bSJed Brown 
49c4762a1bSJed Brown /* Sizes and indices */
50c4762a1bSJed Brown const PetscInt nbus    = 9;         /* Number of network buses */
51c4762a1bSJed Brown const PetscInt ngen    = 3;         /* Number of generators */
52c4762a1bSJed Brown const PetscInt nload   = 3;         /* Number of loads */
53c4762a1bSJed Brown const PetscInt gbus[3] = {0, 1, 2}; /* Buses at which generators are incident */
54c4762a1bSJed Brown const PetscInt lbus[3] = {4, 5, 7}; /* Buses at which loads are incident */
55c4762a1bSJed Brown 
56c4762a1bSJed Brown /* Generator real and reactive powers (found via loadflow) */
57c4762a1bSJed Brown const PetscScalar PG[3] = {0.716786142395021, 1.630000000000000, 0.850000000000000};
58c4762a1bSJed Brown const PetscScalar QG[3] = {0.270702180178785, 0.066120127797275, -0.108402221791588};
59c4762a1bSJed Brown /* Generator constants */
60c4762a1bSJed Brown const PetscScalar H[3]    = {23.64, 6.4, 3.01};       /* Inertia constant */
61c4762a1bSJed Brown const PetscScalar Rs[3]   = {0.0, 0.0, 0.0};          /* Stator Resistance */
62c4762a1bSJed Brown const PetscScalar Xd[3]   = {0.146, 0.8958, 1.3125};  /* d-axis reactance */
63c4762a1bSJed Brown const PetscScalar Xdp[3]  = {0.0608, 0.1198, 0.1813}; /* d-axis transient reactance */
64c4762a1bSJed Brown 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 */
65c4762a1bSJed Brown const PetscScalar Xqp[3]  = {0.0969, 0.1969, 0.25};   /* q-axis transient reactance */
66c4762a1bSJed Brown const PetscScalar Td0p[3] = {8.96, 6.0, 5.89};        /* d-axis open circuit time constant */
67c4762a1bSJed Brown const PetscScalar Tq0p[3] = {0.31, 0.535, 0.6};       /* q-axis open circuit time constant */
68c4762a1bSJed Brown PetscScalar       M[3];                               /* M = 2*H/w_s */
69c4762a1bSJed Brown PetscScalar       D[3];                               /* D = 0.1*M */
70c4762a1bSJed Brown 
71c4762a1bSJed Brown PetscScalar TM[3]; /* Mechanical Torque */
72c4762a1bSJed Brown /* Exciter system constants */
73c4762a1bSJed Brown const PetscScalar KA[3]    = {20.0, 20.0, 20.0};    /* Voltage regulartor gain constant */
74c4762a1bSJed Brown const PetscScalar TA[3]    = {0.2, 0.2, 0.2};       /* Voltage regulator time constant */
75c4762a1bSJed Brown const PetscScalar KE[3]    = {1.0, 1.0, 1.0};       /* Exciter gain constant */
76c4762a1bSJed Brown const PetscScalar TE[3]    = {0.314, 0.314, 0.314}; /* Exciter time constant */
77c4762a1bSJed Brown const PetscScalar KF[3]    = {0.063, 0.063, 0.063}; /* Feedback stabilizer gain constant */
78c4762a1bSJed Brown const PetscScalar TF[3]    = {0.35, 0.35, 0.35};    /* Feedback stabilizer time constant */
79c4762a1bSJed Brown const PetscScalar k1[3]    = {0.0039, 0.0039, 0.0039};
80c4762a1bSJed Brown const PetscScalar k2[3]    = {1.555, 1.555, 1.555}; /* k1 and k2 for calculating the saturation function SE = k1*exp(k2*Efd) */
81c4762a1bSJed Brown const PetscScalar VRMIN[3] = {-4.0, -4.0, -4.0};
82c4762a1bSJed Brown const PetscScalar VRMAX[3] = {7.0, 7.0, 7.0};
83c4762a1bSJed Brown PetscInt          VRatmin[3];
84c4762a1bSJed Brown PetscInt          VRatmax[3];
85c4762a1bSJed Brown 
86c4762a1bSJed Brown PetscScalar Vref[3];
87c4762a1bSJed Brown /* Load constants
88c4762a1bSJed Brown   We use a composite load model that describes the load and reactive powers at each time instant as follows
89c4762a1bSJed Brown   P(t) = \sum\limits_{i=0}^ld_nsegsp \ld_alphap_i*P_D0(\frac{V_m(t)}{V_m0})^\ld_betap_i
90c4762a1bSJed Brown   Q(t) = \sum\limits_{i=0}^ld_nsegsq \ld_alphaq_i*Q_D0(\frac{V_m(t)}{V_m0})^\ld_betaq_i
91c4762a1bSJed Brown   where
92c4762a1bSJed Brown     ld_nsegsp,ld_nsegsq - Number of individual load models for real and reactive power loads
93c4762a1bSJed Brown     ld_alphap,ld_alphap - Percentage contribution (weights) or loads
94c4762a1bSJed Brown     P_D0                - Real power load
95c4762a1bSJed Brown     Q_D0                - Reactive power load
96c4762a1bSJed Brown     V_m(t)              - Voltage magnitude at time t
97c4762a1bSJed Brown     V_m0                - Voltage magnitude at t = 0
98c4762a1bSJed Brown     ld_betap, ld_betaq  - exponents describing the load model for real and reactive part
99c4762a1bSJed Brown 
100c4762a1bSJed Brown     Note: All loads have the same characteristic currently.
101c4762a1bSJed Brown */
102c4762a1bSJed Brown const PetscScalar PD0[3]       = {1.25, 0.9, 1.0};
103c4762a1bSJed Brown const PetscScalar QD0[3]       = {0.5, 0.3, 0.35};
104c4762a1bSJed Brown const PetscInt    ld_nsegsp[3] = {3, 3, 3};
105c4762a1bSJed Brown const PetscScalar ld_alphap[3] = {1.0, 0.0, 0.0};
106c4762a1bSJed Brown const PetscScalar ld_betap[3]  = {2.0, 1.0, 0.0};
107c4762a1bSJed Brown const PetscInt    ld_nsegsq[3] = {3, 3, 3};
108c4762a1bSJed Brown const PetscScalar ld_alphaq[3] = {1.0, 0.0, 0.0};
109c4762a1bSJed Brown const PetscScalar ld_betaq[3]  = {2.0, 1.0, 0.0};
110c4762a1bSJed Brown 
111c4762a1bSJed Brown typedef struct {
112c4762a1bSJed Brown   DM          dmgen, dmnet;        /* DMs to manage generator and network subsystem */
113c4762a1bSJed Brown   DM          dmpgrid;             /* Composite DM to manage the entire power grid */
114c4762a1bSJed Brown   Mat         Ybus;                /* Network admittance matrix */
115c4762a1bSJed Brown   Vec         V0;                  /* Initial voltage vector (Power flow solution) */
116c4762a1bSJed Brown   PetscReal   tfaulton, tfaultoff; /* Fault on and off times */
117c4762a1bSJed Brown   PetscInt    faultbus;            /* Fault bus */
118c4762a1bSJed Brown   PetscScalar Rfault;
119c4762a1bSJed Brown   PetscReal   t0, tmax;
120c4762a1bSJed Brown   PetscInt    neqs_gen, neqs_net, neqs_pgrid;
121c4762a1bSJed Brown   Mat         Sol; /* Matrix to save solution at each time step */
122c4762a1bSJed Brown   PetscInt    stepnum;
123c4762a1bSJed Brown   PetscReal   t;
124c4762a1bSJed Brown   SNES        snes_alg;
125c4762a1bSJed Brown   IS          is_diff;      /* indices for differential equations */
126c4762a1bSJed Brown   IS          is_alg;       /* indices for algebraic equations */
127c4762a1bSJed Brown   PetscBool   setisdiff;    /* TS computes truncation error based only on the differential variables */
128c4762a1bSJed Brown   PetscBool   semiexplicit; /* If the flag is set then a semi-explicit method is used using TSRK */
129c4762a1bSJed Brown } Userctx;
130c4762a1bSJed Brown 
131c4762a1bSJed Brown /*
132c4762a1bSJed Brown   The first two events are for fault on and off, respectively. The following events are
133c4762a1bSJed Brown   to check the min/max limits on the state variable VR. A non windup limiter is used for
134c4762a1bSJed Brown   the VR limits.
135c4762a1bSJed Brown */
EventFunction(TS ts,PetscReal t,Vec X,PetscReal * fvalue,PetscCtx ctx)136*2a8381b2SBarry Smith PetscErrorCode EventFunction(TS ts, PetscReal t, Vec X, PetscReal *fvalue, PetscCtx ctx)
137d71ae5a4SJacob Faibussowitsch {
138c4762a1bSJed Brown   Userctx           *user = (Userctx *)ctx;
139c4762a1bSJed Brown   Vec                Xgen, Xnet;
140c4762a1bSJed Brown   PetscInt           i, idx = 0;
141c4762a1bSJed Brown   const PetscScalar *xgen, *xnet;
142c4762a1bSJed Brown   PetscScalar        Efd, RF, VR, Vr, Vi, Vm;
143c4762a1bSJed Brown 
144c4762a1bSJed Brown   PetscFunctionBegin;
1459566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetLocalVectors(user->dmpgrid, &Xgen, &Xnet));
1469566063dSJacob Faibussowitsch   PetscCall(DMCompositeScatter(user->dmpgrid, X, Xgen, Xnet));
147c4762a1bSJed Brown 
1489566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xgen, &xgen));
1499566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xnet, &xnet));
150c4762a1bSJed Brown 
151c4762a1bSJed Brown   /* Event for fault-on time */
152c4762a1bSJed Brown   fvalue[0] = t - user->tfaulton;
153c4762a1bSJed Brown   /* Event for fault-off time */
154c4762a1bSJed Brown   fvalue[1] = t - user->tfaultoff;
155c4762a1bSJed Brown 
156c4762a1bSJed Brown   for (i = 0; i < ngen; i++) {
157c4762a1bSJed Brown     Efd = xgen[idx + 6];
158c4762a1bSJed Brown     RF  = xgen[idx + 7];
159c4762a1bSJed Brown     VR  = xgen[idx + 8];
160c4762a1bSJed Brown 
161c4762a1bSJed Brown     Vr = xnet[2 * gbus[i]];     /* Real part of generator terminal voltage */
162c4762a1bSJed Brown     Vi = xnet[2 * gbus[i] + 1]; /* Imaginary part of the generator terminal voltage */
163c4762a1bSJed Brown     Vm = PetscSqrtScalar(Vr * Vr + Vi * Vi);
164c4762a1bSJed Brown 
165c4762a1bSJed Brown     if (!VRatmax[i]) {
166fa9584fbSIlya Fursov       fvalue[2 + 2 * i] = PetscRealPart(VRMAX[i] - VR);
167c4762a1bSJed Brown     } else {
168fa9584fbSIlya Fursov       fvalue[2 + 2 * i] = PetscRealPart((VR - KA[i] * RF + KA[i] * KF[i] * Efd / TF[i] - KA[i] * (Vref[i] - Vm)) / TA[i]);
169c4762a1bSJed Brown     }
170c4762a1bSJed Brown     if (!VRatmin[i]) {
171fa9584fbSIlya Fursov       fvalue[2 + 2 * i + 1] = PetscRealPart(VRMIN[i] - VR);
172c4762a1bSJed Brown     } else {
173fa9584fbSIlya Fursov       fvalue[2 + 2 * i + 1] = PetscRealPart((VR - KA[i] * RF + KA[i] * KF[i] * Efd / TF[i] - KA[i] * (Vref[i] - Vm)) / TA[i]);
174c4762a1bSJed Brown     }
175c4762a1bSJed Brown     idx = idx + 9;
176c4762a1bSJed Brown   }
1779566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xgen, &xgen));
1789566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xnet, &xnet));
179c4762a1bSJed Brown 
1809566063dSJacob Faibussowitsch   PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid, &Xgen, &Xnet));
1813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
182c4762a1bSJed Brown }
183c4762a1bSJed Brown 
PostEventFunction(TS ts,PetscInt nevents,PetscInt event_list[],PetscReal t,Vec X,PetscBool forwardsolve,PetscCtx ctx)184*2a8381b2SBarry Smith PetscErrorCode PostEventFunction(TS ts, PetscInt nevents, PetscInt event_list[], PetscReal t, Vec X, PetscBool forwardsolve, PetscCtx ctx)
185d71ae5a4SJacob Faibussowitsch {
186c4762a1bSJed Brown   Userctx     *user = (Userctx *)ctx;
187c4762a1bSJed Brown   Vec          Xgen, Xnet;
188c4762a1bSJed Brown   PetscScalar *xgen, *xnet;
189c4762a1bSJed Brown   PetscInt     row_loc, col_loc;
190c4762a1bSJed Brown   PetscScalar  val;
191c4762a1bSJed Brown   PetscInt     i, idx = 0, event_num;
192c4762a1bSJed Brown   PetscScalar  fvalue;
193c4762a1bSJed Brown   PetscScalar  Efd, RF, VR;
194c4762a1bSJed Brown   PetscScalar  Vr, Vi, Vm;
195c4762a1bSJed Brown 
196c4762a1bSJed Brown   PetscFunctionBegin;
1979566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetLocalVectors(user->dmpgrid, &Xgen, &Xnet));
1989566063dSJacob Faibussowitsch   PetscCall(DMCompositeScatter(user->dmpgrid, X, Xgen, Xnet));
199c4762a1bSJed Brown 
2009566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Xgen, &xgen));
2019566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Xnet, &xnet));
202c4762a1bSJed Brown 
203c4762a1bSJed Brown   for (i = 0; i < nevents; i++) {
204c4762a1bSJed Brown     if (event_list[i] == 0) {
205c4762a1bSJed Brown       /* Apply disturbance - resistive fault at user->faultbus */
206c4762a1bSJed Brown       /* This is done by adding shunt conductance to the diagonal location
207c4762a1bSJed Brown          in the Ybus matrix */
2089371c9d4SSatish Balay       row_loc = 2 * user->faultbus;
2099371c9d4SSatish Balay       col_loc = 2 * user->faultbus + 1; /* Location for G */
210c4762a1bSJed Brown       val     = 1 / user->Rfault;
2119566063dSJacob Faibussowitsch       PetscCall(MatSetValues(user->Ybus, 1, &row_loc, 1, &col_loc, &val, ADD_VALUES));
2129371c9d4SSatish Balay       row_loc = 2 * user->faultbus + 1;
2139371c9d4SSatish Balay       col_loc = 2 * user->faultbus; /* Location for G */
214c4762a1bSJed Brown       val     = 1 / user->Rfault;
2159566063dSJacob Faibussowitsch       PetscCall(MatSetValues(user->Ybus, 1, &row_loc, 1, &col_loc, &val, ADD_VALUES));
216c4762a1bSJed Brown 
2179566063dSJacob Faibussowitsch       PetscCall(MatAssemblyBegin(user->Ybus, MAT_FINAL_ASSEMBLY));
2189566063dSJacob Faibussowitsch       PetscCall(MatAssemblyEnd(user->Ybus, MAT_FINAL_ASSEMBLY));
219c4762a1bSJed Brown 
220c4762a1bSJed Brown       /* Solve the algebraic equations */
2219566063dSJacob Faibussowitsch       PetscCall(SNESSolve(user->snes_alg, NULL, X));
222c4762a1bSJed Brown     } else if (event_list[i] == 1) {
223c4762a1bSJed Brown       /* Remove the fault */
2249371c9d4SSatish Balay       row_loc = 2 * user->faultbus;
2259371c9d4SSatish Balay       col_loc = 2 * user->faultbus + 1;
226c4762a1bSJed Brown       val     = -1 / user->Rfault;
2279566063dSJacob Faibussowitsch       PetscCall(MatSetValues(user->Ybus, 1, &row_loc, 1, &col_loc, &val, ADD_VALUES));
2289371c9d4SSatish Balay       row_loc = 2 * user->faultbus + 1;
2299371c9d4SSatish Balay       col_loc = 2 * user->faultbus;
230c4762a1bSJed Brown       val     = -1 / user->Rfault;
2319566063dSJacob Faibussowitsch       PetscCall(MatSetValues(user->Ybus, 1, &row_loc, 1, &col_loc, &val, ADD_VALUES));
232c4762a1bSJed Brown 
2339566063dSJacob Faibussowitsch       PetscCall(MatAssemblyBegin(user->Ybus, MAT_FINAL_ASSEMBLY));
2349566063dSJacob Faibussowitsch       PetscCall(MatAssemblyEnd(user->Ybus, MAT_FINAL_ASSEMBLY));
235c4762a1bSJed Brown 
236c4762a1bSJed Brown       /* Solve the algebraic equations */
2379566063dSJacob Faibussowitsch       PetscCall(SNESSolve(user->snes_alg, NULL, X));
238c4762a1bSJed Brown 
239c4762a1bSJed Brown       /* Check the VR derivatives and reset flags if needed */
240c4762a1bSJed Brown       for (i = 0; i < ngen; i++) {
241c4762a1bSJed Brown         Efd = xgen[idx + 6];
242c4762a1bSJed Brown         RF  = xgen[idx + 7];
243c4762a1bSJed Brown         VR  = xgen[idx + 8];
244c4762a1bSJed Brown 
245c4762a1bSJed Brown         Vr = xnet[2 * gbus[i]];     /* Real part of generator terminal voltage */
246c4762a1bSJed Brown         Vi = xnet[2 * gbus[i] + 1]; /* Imaginary part of the generator terminal voltage */
247c4762a1bSJed Brown         Vm = PetscSqrtScalar(Vr * Vr + Vi * Vi);
248c4762a1bSJed Brown 
249c4762a1bSJed Brown         if (VRatmax[i]) {
250c4762a1bSJed Brown           fvalue = (VR - KA[i] * RF + KA[i] * KF[i] * Efd / TF[i] - KA[i] * (Vref[i] - Vm)) / TA[i];
251c4762a1bSJed Brown           if (fvalue < 0) {
252c4762a1bSJed Brown             VRatmax[i] = 0;
25363a3b9bcSJacob Faibussowitsch             PetscCall(PetscPrintf(PETSC_COMM_SELF, "VR[%" PetscInt_FMT "]: dVR_dt went negative on fault clearing at time %g\n", i, (double)t));
254c4762a1bSJed Brown           }
255c4762a1bSJed Brown         }
256c4762a1bSJed Brown         if (VRatmin[i]) {
257c4762a1bSJed Brown           fvalue = (VR - KA[i] * RF + KA[i] * KF[i] * Efd / TF[i] - KA[i] * (Vref[i] - Vm)) / TA[i];
258c4762a1bSJed Brown 
259c4762a1bSJed Brown           if (fvalue > 0) {
260c4762a1bSJed Brown             VRatmin[i] = 0;
26163a3b9bcSJacob Faibussowitsch             PetscCall(PetscPrintf(PETSC_COMM_SELF, "VR[%" PetscInt_FMT "]: dVR_dt went positive on fault clearing at time %g\n", i, (double)t));
262c4762a1bSJed Brown           }
263c4762a1bSJed Brown         }
264c4762a1bSJed Brown         idx = idx + 9;
265c4762a1bSJed Brown       }
266c4762a1bSJed Brown     } else {
267c4762a1bSJed Brown       idx       = (event_list[i] - 2) / 2;
268c4762a1bSJed Brown       event_num = (event_list[i] - 2) % 2;
269c4762a1bSJed Brown       if (event_num == 0) { /* Max VR */
270c4762a1bSJed Brown         if (!VRatmax[idx]) {
271c4762a1bSJed Brown           VRatmax[idx] = 1;
27263a3b9bcSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "VR[%" PetscInt_FMT "]: hit upper limit at time %g\n", idx, (double)t));
2739371c9d4SSatish Balay         } else {
274c4762a1bSJed Brown           VRatmax[idx] = 0;
27563a3b9bcSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "VR[%" PetscInt_FMT "]: freeing variable as dVR_dt is negative at time %g\n", idx, (double)t));
276c4762a1bSJed Brown         }
277c4762a1bSJed Brown       } else {
278c4762a1bSJed Brown         if (!VRatmin[idx]) {
279c4762a1bSJed Brown           VRatmin[idx] = 1;
28063a3b9bcSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "VR[%" PetscInt_FMT "]: hit lower limit at time %g\n", idx, (double)t));
2819371c9d4SSatish Balay         } else {
282c4762a1bSJed Brown           VRatmin[idx] = 0;
28363a3b9bcSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "VR[%" PetscInt_FMT "]: freeing variable as dVR_dt is positive at time %g\n", idx, (double)t));
284c4762a1bSJed Brown         }
285c4762a1bSJed Brown       }
286c4762a1bSJed Brown     }
287c4762a1bSJed Brown   }
288c4762a1bSJed Brown 
2899566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Xgen, &xgen));
2909566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Xnet, &xnet));
291c4762a1bSJed Brown 
2929566063dSJacob Faibussowitsch   PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid, &Xgen, &Xnet));
2933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
294c4762a1bSJed Brown }
295c4762a1bSJed Brown 
296c4762a1bSJed Brown /* Converts from machine frame (dq) to network (phase a real,imag) reference frame */
dq2ri(PetscScalar Fd,PetscScalar Fq,PetscScalar delta,PetscScalar * Fr,PetscScalar * Fi)297d71ae5a4SJacob Faibussowitsch PetscErrorCode dq2ri(PetscScalar Fd, PetscScalar Fq, PetscScalar delta, PetscScalar *Fr, PetscScalar *Fi)
298d71ae5a4SJacob Faibussowitsch {
299c4762a1bSJed Brown   PetscFunctionBegin;
300c4762a1bSJed Brown   *Fr = Fd * PetscSinScalar(delta) + Fq * PetscCosScalar(delta);
301c4762a1bSJed Brown   *Fi = -Fd * PetscCosScalar(delta) + Fq * PetscSinScalar(delta);
3023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
303c4762a1bSJed Brown }
304c4762a1bSJed Brown 
305c4762a1bSJed Brown /* Converts from network frame ([phase a real,imag) to machine (dq) reference frame */
ri2dq(PetscScalar Fr,PetscScalar Fi,PetscScalar delta,PetscScalar * Fd,PetscScalar * Fq)306d71ae5a4SJacob Faibussowitsch PetscErrorCode ri2dq(PetscScalar Fr, PetscScalar Fi, PetscScalar delta, PetscScalar *Fd, PetscScalar *Fq)
307d71ae5a4SJacob Faibussowitsch {
308c4762a1bSJed Brown   PetscFunctionBegin;
309c4762a1bSJed Brown   *Fd = Fr * PetscSinScalar(delta) - Fi * PetscCosScalar(delta);
310c4762a1bSJed Brown   *Fq = Fr * PetscCosScalar(delta) + Fi * PetscSinScalar(delta);
3113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
312c4762a1bSJed Brown }
313c4762a1bSJed Brown 
314c4762a1bSJed Brown /* Saves the solution at each time to a matrix */
SaveSolution(TS ts)315d71ae5a4SJacob Faibussowitsch PetscErrorCode SaveSolution(TS ts)
316d71ae5a4SJacob Faibussowitsch {
317c4762a1bSJed Brown   Userctx           *user;
318c4762a1bSJed Brown   Vec                X;
319c4762a1bSJed Brown   const PetscScalar *x;
320c4762a1bSJed Brown   PetscScalar       *mat;
321c4762a1bSJed Brown   PetscInt           idx;
322c4762a1bSJed Brown   PetscReal          t;
323c4762a1bSJed Brown 
324c4762a1bSJed Brown   PetscFunctionBegin;
3259566063dSJacob Faibussowitsch   PetscCall(TSGetApplicationContext(ts, &user));
3269566063dSJacob Faibussowitsch   PetscCall(TSGetTime(ts, &t));
3279566063dSJacob Faibussowitsch   PetscCall(TSGetSolution(ts, &X));
328c4762a1bSJed Brown   idx = user->stepnum * (user->neqs_pgrid + 1);
3299566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArray(user->Sol, &mat));
3309566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
331c4762a1bSJed Brown   mat[idx] = t;
3329566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(mat + idx + 1, x, user->neqs_pgrid));
3339566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArray(user->Sol, &mat));
3349566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
335c4762a1bSJed Brown   user->stepnum++;
3363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
337c4762a1bSJed Brown }
338c4762a1bSJed Brown 
SetInitialGuess(Vec X,Userctx * user)339d71ae5a4SJacob Faibussowitsch PetscErrorCode SetInitialGuess(Vec X, Userctx *user)
340d71ae5a4SJacob Faibussowitsch {
341c4762a1bSJed Brown   Vec                Xgen, Xnet;
34205e808d2SBarry Smith   PetscScalar       *xgen;
34305e808d2SBarry Smith   const PetscScalar *xnet;
344c4762a1bSJed Brown   PetscInt           i, idx = 0;
345c4762a1bSJed Brown   PetscScalar        Vr, Vi, IGr, IGi, Vm, Vm2;
346c4762a1bSJed Brown   PetscScalar        Eqp, Edp, delta;
347c4762a1bSJed Brown   PetscScalar        Efd, RF, VR; /* Exciter variables */
348c4762a1bSJed Brown   PetscScalar        Id, Iq;      /* Generator dq axis currents */
349c4762a1bSJed Brown   PetscScalar        theta, Vd, Vq, SE;
350c4762a1bSJed Brown 
351c4762a1bSJed Brown   PetscFunctionBegin;
3529371c9d4SSatish Balay   M[0] = 2 * H[0] / w_s;
3539371c9d4SSatish Balay   M[1] = 2 * H[1] / w_s;
3549371c9d4SSatish Balay   M[2] = 2 * H[2] / w_s;
3559371c9d4SSatish Balay   D[0] = 0.1 * M[0];
3569371c9d4SSatish Balay   D[1] = 0.1 * M[1];
3579371c9d4SSatish Balay   D[2] = 0.1 * M[2];
358c4762a1bSJed Brown 
3599566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetLocalVectors(user->dmpgrid, &Xgen, &Xnet));
360c4762a1bSJed Brown 
361c4762a1bSJed Brown   /* Network subsystem initialization */
3629566063dSJacob Faibussowitsch   PetscCall(VecCopy(user->V0, Xnet));
363c4762a1bSJed Brown 
364c4762a1bSJed Brown   /* Generator subsystem initialization */
3659566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(Xgen, &xgen));
3669566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xnet, &xnet));
367c4762a1bSJed Brown 
368c4762a1bSJed Brown   for (i = 0; i < ngen; i++) {
369c4762a1bSJed Brown     Vr  = xnet[2 * gbus[i]];     /* Real part of generator terminal voltage */
370c4762a1bSJed Brown     Vi  = xnet[2 * gbus[i] + 1]; /* Imaginary part of the generator terminal voltage */
3719371c9d4SSatish Balay     Vm  = PetscSqrtScalar(Vr * Vr + Vi * Vi);
3729371c9d4SSatish Balay     Vm2 = Vm * Vm;
373c4762a1bSJed Brown     IGr = (Vr * PG[i] + Vi * QG[i]) / Vm2;
374c4762a1bSJed Brown     IGi = (Vi * PG[i] - Vr * QG[i]) / Vm2;
375c4762a1bSJed Brown 
376c4762a1bSJed Brown     delta = PetscAtan2Real(Vi + Xq[i] * IGr, Vr - Xq[i] * IGi); /* Machine angle */
377c4762a1bSJed Brown 
378c4762a1bSJed Brown     theta = PETSC_PI / 2.0 - delta;
379c4762a1bSJed Brown 
380c4762a1bSJed Brown     Id = IGr * PetscCosScalar(theta) - IGi * PetscSinScalar(theta); /* d-axis stator current */
381c4762a1bSJed Brown     Iq = IGr * PetscSinScalar(theta) + IGi * PetscCosScalar(theta); /* q-axis stator current */
382c4762a1bSJed Brown 
383c4762a1bSJed Brown     Vd = Vr * PetscCosScalar(theta) - Vi * PetscSinScalar(theta);
384c4762a1bSJed Brown     Vq = Vr * PetscSinScalar(theta) + Vi * PetscCosScalar(theta);
385c4762a1bSJed Brown 
386c4762a1bSJed Brown     Edp = Vd + Rs[i] * Id - Xqp[i] * Iq; /* d-axis transient EMF */
387c4762a1bSJed Brown     Eqp = Vq + Rs[i] * Iq + Xdp[i] * Id; /* q-axis transient EMF */
388c4762a1bSJed Brown 
389c4762a1bSJed Brown     TM[i] = PG[i];
390c4762a1bSJed Brown 
391c4762a1bSJed Brown     /* The generator variables are ordered as [Eqp,Edp,delta,w,Id,Iq] */
392c4762a1bSJed Brown     xgen[idx]     = Eqp;
393c4762a1bSJed Brown     xgen[idx + 1] = Edp;
394c4762a1bSJed Brown     xgen[idx + 2] = delta;
395c4762a1bSJed Brown     xgen[idx + 3] = w_s;
396c4762a1bSJed Brown 
397c4762a1bSJed Brown     idx = idx + 4;
398c4762a1bSJed Brown 
399c4762a1bSJed Brown     xgen[idx]     = Id;
400c4762a1bSJed Brown     xgen[idx + 1] = Iq;
401c4762a1bSJed Brown 
402c4762a1bSJed Brown     idx = idx + 2;
403c4762a1bSJed Brown 
404c4762a1bSJed Brown     /* Exciter */
405c4762a1bSJed Brown     Efd = Eqp + (Xd[i] - Xdp[i]) * Id;
406c4762a1bSJed Brown     SE  = k1[i] * PetscExpScalar(k2[i] * Efd);
407c4762a1bSJed Brown     VR  = KE[i] * Efd + SE;
408c4762a1bSJed Brown     RF  = KF[i] * Efd / TF[i];
409c4762a1bSJed Brown 
410c4762a1bSJed Brown     xgen[idx]     = Efd;
411c4762a1bSJed Brown     xgen[idx + 1] = RF;
412c4762a1bSJed Brown     xgen[idx + 2] = VR;
413c4762a1bSJed Brown 
414c4762a1bSJed Brown     Vref[i] = Vm + (VR / KA[i]);
415c4762a1bSJed Brown 
416c4762a1bSJed Brown     VRatmin[i] = VRatmax[i] = 0;
417c4762a1bSJed Brown 
418c4762a1bSJed Brown     idx = idx + 3;
419c4762a1bSJed Brown   }
4209566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(Xgen, &xgen));
4219566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xnet, &xnet));
422c4762a1bSJed Brown 
4239566063dSJacob Faibussowitsch   /* PetscCall(VecView(Xgen,0)); */
4249566063dSJacob Faibussowitsch   PetscCall(DMCompositeGather(user->dmpgrid, INSERT_VALUES, X, Xgen, Xnet));
4259566063dSJacob Faibussowitsch   PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid, &Xgen, &Xnet));
4263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
427c4762a1bSJed Brown }
428c4762a1bSJed Brown 
429c4762a1bSJed Brown /* Computes F = [f(x,y);g(x,y)] */
ResidualFunction(Vec X,Vec F,Userctx * user)430d71ae5a4SJacob Faibussowitsch PetscErrorCode ResidualFunction(Vec X, Vec F, Userctx *user)
431d71ae5a4SJacob Faibussowitsch {
432c4762a1bSJed Brown   Vec                Xgen, Xnet, Fgen, Fnet;
43305e808d2SBarry Smith   const PetscScalar *xgen, *xnet;
43405e808d2SBarry Smith   PetscScalar       *fgen, *fnet;
435c4762a1bSJed Brown   PetscInt           i, idx = 0;
436c4762a1bSJed Brown   PetscScalar        Vr, Vi, Vm, Vm2;
437c4762a1bSJed Brown   PetscScalar        Eqp, Edp, delta, w; /* Generator variables */
438c4762a1bSJed Brown   PetscScalar        Efd, RF, VR;        /* Exciter variables */
439c4762a1bSJed Brown   PetscScalar        Id, Iq;             /* Generator dq axis currents */
440c4762a1bSJed Brown   PetscScalar        Vd, Vq, SE;
441c4762a1bSJed Brown   PetscScalar        IGr, IGi, IDr, IDi;
442c4762a1bSJed Brown   PetscScalar        Zdq_inv[4], det;
443c4762a1bSJed Brown   PetscScalar        PD, QD, Vm0, *v0;
444c4762a1bSJed Brown   PetscInt           k;
445c4762a1bSJed Brown 
446c4762a1bSJed Brown   PetscFunctionBegin;
4479566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(F));
4489566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetLocalVectors(user->dmpgrid, &Xgen, &Xnet));
4499566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetLocalVectors(user->dmpgrid, &Fgen, &Fnet));
4509566063dSJacob Faibussowitsch   PetscCall(DMCompositeScatter(user->dmpgrid, X, Xgen, Xnet));
4519566063dSJacob Faibussowitsch   PetscCall(DMCompositeScatter(user->dmpgrid, F, Fgen, Fnet));
452c4762a1bSJed Brown 
453c4762a1bSJed Brown   /* Network current balance residual IG + Y*V + IL = 0. Only YV is added here.
454c4762a1bSJed Brown      The generator current injection, IG, and load current injection, ID are added later
455c4762a1bSJed Brown   */
456c4762a1bSJed Brown   /* Note that the values in Ybus are stored assuming the imaginary current balance
457c4762a1bSJed Brown      equation is ordered first followed by real current balance equation for each bus.
458c4762a1bSJed Brown      Thus imaginary current contribution goes in location 2*i, and
459c4762a1bSJed Brown      real current contribution in 2*i+1
460c4762a1bSJed Brown   */
4619566063dSJacob Faibussowitsch   PetscCall(MatMult(user->Ybus, Xnet, Fnet));
462c4762a1bSJed Brown 
4639566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xgen, &xgen));
4649566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xnet, &xnet));
4659566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(Fgen, &fgen));
4669566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(Fnet, &fnet));
467c4762a1bSJed Brown 
468c4762a1bSJed Brown   /* Generator subsystem */
469c4762a1bSJed Brown   for (i = 0; i < ngen; i++) {
470c4762a1bSJed Brown     Eqp   = xgen[idx];
471c4762a1bSJed Brown     Edp   = xgen[idx + 1];
472c4762a1bSJed Brown     delta = xgen[idx + 2];
473c4762a1bSJed Brown     w     = xgen[idx + 3];
474c4762a1bSJed Brown     Id    = xgen[idx + 4];
475c4762a1bSJed Brown     Iq    = xgen[idx + 5];
476c4762a1bSJed Brown     Efd   = xgen[idx + 6];
477c4762a1bSJed Brown     RF    = xgen[idx + 7];
478c4762a1bSJed Brown     VR    = xgen[idx + 8];
479c4762a1bSJed Brown 
480c4762a1bSJed Brown     /* Generator differential equations */
481c4762a1bSJed Brown     fgen[idx]     = (-Eqp - (Xd[i] - Xdp[i]) * Id + Efd) / Td0p[i];
482c4762a1bSJed Brown     fgen[idx + 1] = (-Edp + (Xq[i] - Xqp[i]) * Iq) / Tq0p[i];
483c4762a1bSJed Brown     fgen[idx + 2] = w - w_s;
484c4762a1bSJed Brown     fgen[idx + 3] = (TM[i] - Edp * Id - Eqp * Iq - (Xqp[i] - Xdp[i]) * Id * Iq - D[i] * (w - w_s)) / M[i];
485c4762a1bSJed Brown 
486c4762a1bSJed Brown     Vr = xnet[2 * gbus[i]];     /* Real part of generator terminal voltage */
487c4762a1bSJed Brown     Vi = xnet[2 * gbus[i] + 1]; /* Imaginary part of the generator terminal voltage */
488c4762a1bSJed Brown 
4899566063dSJacob Faibussowitsch     PetscCall(ri2dq(Vr, Vi, delta, &Vd, &Vq));
490c4762a1bSJed Brown     /* Algebraic equations for stator currents */
491c4762a1bSJed Brown     det = Rs[i] * Rs[i] + Xdp[i] * Xqp[i];
492c4762a1bSJed Brown 
493c4762a1bSJed Brown     Zdq_inv[0] = Rs[i] / det;
494c4762a1bSJed Brown     Zdq_inv[1] = Xqp[i] / det;
495c4762a1bSJed Brown     Zdq_inv[2] = -Xdp[i] / det;
496c4762a1bSJed Brown     Zdq_inv[3] = Rs[i] / det;
497c4762a1bSJed Brown 
498c4762a1bSJed Brown     fgen[idx + 4] = Zdq_inv[0] * (-Edp + Vd) + Zdq_inv[1] * (-Eqp + Vq) + Id;
499c4762a1bSJed Brown     fgen[idx + 5] = Zdq_inv[2] * (-Edp + Vd) + Zdq_inv[3] * (-Eqp + Vq) + Iq;
500c4762a1bSJed Brown 
501c4762a1bSJed Brown     /* Add generator current injection to network */
5029566063dSJacob Faibussowitsch     PetscCall(dq2ri(Id, Iq, delta, &IGr, &IGi));
503c4762a1bSJed Brown 
504c4762a1bSJed Brown     fnet[2 * gbus[i]] -= IGi;
505c4762a1bSJed Brown     fnet[2 * gbus[i] + 1] -= IGr;
506c4762a1bSJed Brown 
507c4762a1bSJed Brown     Vm = PetscSqrtScalar(Vd * Vd + Vq * Vq);
508c4762a1bSJed Brown 
509c4762a1bSJed Brown     SE = k1[i] * PetscExpScalar(k2[i] * Efd);
510c4762a1bSJed Brown 
511c4762a1bSJed Brown     /* Exciter differential equations */
512c4762a1bSJed Brown     fgen[idx + 6] = (-KE[i] * Efd - SE + VR) / TE[i];
513c4762a1bSJed Brown     fgen[idx + 7] = (-RF + KF[i] * Efd / TF[i]) / TF[i];
514c4762a1bSJed Brown     if (VRatmax[i]) fgen[idx + 8] = VR - VRMAX[i];
515c4762a1bSJed Brown     else if (VRatmin[i]) fgen[idx + 8] = VRMIN[i] - VR;
516c4762a1bSJed Brown     else fgen[idx + 8] = (-VR + KA[i] * RF - KA[i] * KF[i] * Efd / TF[i] + KA[i] * (Vref[i] - Vm)) / TA[i];
517c4762a1bSJed Brown 
518c4762a1bSJed Brown     idx = idx + 9;
519c4762a1bSJed Brown   }
520c4762a1bSJed Brown 
5219566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user->V0, &v0));
522c4762a1bSJed Brown   for (i = 0; i < nload; i++) {
523c4762a1bSJed Brown     Vr  = xnet[2 * lbus[i]];     /* Real part of load bus voltage */
524c4762a1bSJed Brown     Vi  = xnet[2 * lbus[i] + 1]; /* Imaginary part of the load bus voltage */
5259371c9d4SSatish Balay     Vm  = PetscSqrtScalar(Vr * Vr + Vi * Vi);
5269371c9d4SSatish Balay     Vm2 = Vm * Vm;
527c4762a1bSJed Brown     Vm0 = PetscSqrtScalar(v0[2 * lbus[i]] * v0[2 * lbus[i]] + v0[2 * lbus[i] + 1] * v0[2 * lbus[i] + 1]);
528c4762a1bSJed Brown     PD = QD = 0.0;
52957508eceSPierre Jolivet     for (k = 0; k < ld_nsegsp[i]; k++) PD += ld_alphap[k] * PD0[i] * PetscPowScalar(Vm / Vm0, ld_betap[k]);
53057508eceSPierre Jolivet     for (k = 0; k < ld_nsegsq[i]; k++) QD += ld_alphaq[k] * QD0[i] * PetscPowScalar(Vm / Vm0, ld_betaq[k]);
531c4762a1bSJed Brown 
532c4762a1bSJed Brown     /* Load currents */
533c4762a1bSJed Brown     IDr = (PD * Vr + QD * Vi) / Vm2;
534c4762a1bSJed Brown     IDi = (-QD * Vr + PD * Vi) / Vm2;
535c4762a1bSJed Brown 
536c4762a1bSJed Brown     fnet[2 * lbus[i]] += IDi;
537c4762a1bSJed Brown     fnet[2 * lbus[i] + 1] += IDr;
538c4762a1bSJed Brown   }
5399566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user->V0, &v0));
540c4762a1bSJed Brown 
5419566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xgen, &xgen));
5429566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xnet, &xnet));
5439566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(Fgen, &fgen));
5449566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(Fnet, &fnet));
545c4762a1bSJed Brown 
5469566063dSJacob Faibussowitsch   PetscCall(DMCompositeGather(user->dmpgrid, INSERT_VALUES, F, Fgen, Fnet));
5479566063dSJacob Faibussowitsch   PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid, &Xgen, &Xnet));
5489566063dSJacob Faibussowitsch   PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid, &Fgen, &Fnet));
5493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
550c4762a1bSJed Brown }
551c4762a1bSJed Brown 
552c4762a1bSJed Brown /*   f(x,y)
553c4762a1bSJed Brown      g(x,y)
554c4762a1bSJed Brown  */
RHSFunction(TS ts,PetscReal t,Vec X,Vec F,PetscCtx ctx)555*2a8381b2SBarry Smith PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec X, Vec F, PetscCtx ctx)
556d71ae5a4SJacob Faibussowitsch {
557c4762a1bSJed Brown   Userctx *user = (Userctx *)ctx;
558c4762a1bSJed Brown 
559c4762a1bSJed Brown   PetscFunctionBegin;
560c4762a1bSJed Brown   user->t = t;
5619566063dSJacob Faibussowitsch   PetscCall(ResidualFunction(X, F, user));
5623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
563c4762a1bSJed Brown }
564c4762a1bSJed Brown 
565c4762a1bSJed Brown /* f(x,y) - \dot{x}
566c4762a1bSJed Brown      g(x,y) = 0
567c4762a1bSJed Brown  */
IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,PetscCtx ctx)568*2a8381b2SBarry Smith PetscErrorCode IFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, PetscCtx ctx)
569d71ae5a4SJacob Faibussowitsch {
57005e808d2SBarry Smith   PetscScalar       *f;
57105e808d2SBarry Smith   const PetscScalar *xdot;
572c4762a1bSJed Brown   PetscInt           i;
573c4762a1bSJed Brown 
574c4762a1bSJed Brown   PetscFunctionBegin;
5759566063dSJacob Faibussowitsch   PetscCall(RHSFunction(ts, t, X, F, ctx));
5769566063dSJacob Faibussowitsch   PetscCall(VecScale(F, -1.0));
5779566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
5789566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xdot, &xdot));
579c4762a1bSJed Brown   for (i = 0; i < ngen; i++) {
58095a2cb33SBarry Smith     f[9 * i] += xdot[9 * i];
58195a2cb33SBarry Smith     f[9 * i + 1] += xdot[9 * i + 1];
58295a2cb33SBarry Smith     f[9 * i + 2] += xdot[9 * i + 2];
58395a2cb33SBarry Smith     f[9 * i + 3] += xdot[9 * i + 3];
58495a2cb33SBarry Smith     f[9 * i + 6] += xdot[9 * i + 6];
58595a2cb33SBarry Smith     f[9 * i + 7] += xdot[9 * i + 7];
58695a2cb33SBarry Smith     f[9 * i + 8] += xdot[9 * i + 8];
587c4762a1bSJed Brown   }
5889566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
5899566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xdot, &xdot));
5903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
591c4762a1bSJed Brown }
592c4762a1bSJed Brown 
593c4762a1bSJed Brown /* This function is used for solving the algebraic system only during fault on and
594c4762a1bSJed Brown    off times. It computes the entire F and then zeros out the part corresponding to
595c4762a1bSJed Brown    differential equations
596c4762a1bSJed Brown  F = [0;g(y)];
597c4762a1bSJed Brown */
AlgFunction(SNES snes,Vec X,Vec F,PetscCtx ctx)598*2a8381b2SBarry Smith PetscErrorCode AlgFunction(SNES snes, Vec X, Vec F, PetscCtx ctx)
599d71ae5a4SJacob Faibussowitsch {
600c4762a1bSJed Brown   Userctx     *user = (Userctx *)ctx;
601c4762a1bSJed Brown   PetscScalar *f;
602c4762a1bSJed Brown   PetscInt     i;
603c4762a1bSJed Brown 
604c4762a1bSJed Brown   PetscFunctionBegin;
6059566063dSJacob Faibussowitsch   PetscCall(ResidualFunction(X, F, user));
6069566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
607c4762a1bSJed Brown   for (i = 0; i < ngen; i++) {
608c4762a1bSJed Brown     f[9 * i]     = 0;
609c4762a1bSJed Brown     f[9 * i + 1] = 0;
610c4762a1bSJed Brown     f[9 * i + 2] = 0;
611c4762a1bSJed Brown     f[9 * i + 3] = 0;
612c4762a1bSJed Brown     f[9 * i + 6] = 0;
613c4762a1bSJed Brown     f[9 * i + 7] = 0;
614c4762a1bSJed Brown     f[9 * i + 8] = 0;
615c4762a1bSJed Brown   }
6169566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
6173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
618c4762a1bSJed Brown }
619c4762a1bSJed Brown 
PostStage(TS ts,PetscReal t,PetscInt i,Vec * X)620d71ae5a4SJacob Faibussowitsch PetscErrorCode PostStage(TS ts, PetscReal t, PetscInt i, Vec *X)
621d71ae5a4SJacob Faibussowitsch {
622c4762a1bSJed Brown   Userctx *user;
623c4762a1bSJed Brown 
624c4762a1bSJed Brown   PetscFunctionBegin;
6259566063dSJacob Faibussowitsch   PetscCall(TSGetApplicationContext(ts, &user));
6269566063dSJacob Faibussowitsch   PetscCall(SNESSolve(user->snes_alg, NULL, X[i]));
6273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
628c4762a1bSJed Brown }
629c4762a1bSJed Brown 
PostEvaluate(TS ts)630d71ae5a4SJacob Faibussowitsch PetscErrorCode PostEvaluate(TS ts)
631d71ae5a4SJacob Faibussowitsch {
632c4762a1bSJed Brown   Userctx *user;
633c4762a1bSJed Brown   Vec      X;
634c4762a1bSJed Brown 
635c4762a1bSJed Brown   PetscFunctionBegin;
6369566063dSJacob Faibussowitsch   PetscCall(TSGetApplicationContext(ts, &user));
6379566063dSJacob Faibussowitsch   PetscCall(TSGetSolution(ts, &X));
6389566063dSJacob Faibussowitsch   PetscCall(SNESSolve(user->snes_alg, NULL, X));
6393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
640c4762a1bSJed Brown }
641c4762a1bSJed Brown 
PreallocateJacobian(Mat J,Userctx * user)642d71ae5a4SJacob Faibussowitsch PetscErrorCode PreallocateJacobian(Mat J, Userctx *user)
643d71ae5a4SJacob Faibussowitsch {
644c4762a1bSJed Brown   PetscInt *d_nnz;
645c4762a1bSJed Brown   PetscInt  i, idx = 0, start = 0;
646c4762a1bSJed Brown   PetscInt  ncols;
647c4762a1bSJed Brown 
648c4762a1bSJed Brown   PetscFunctionBegin;
6499566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(user->neqs_pgrid, &d_nnz));
650c4762a1bSJed Brown   for (i = 0; i < user->neqs_pgrid; i++) d_nnz[i] = 0;
651c4762a1bSJed Brown   /* Generator subsystem */
652c4762a1bSJed Brown   for (i = 0; i < ngen; i++) {
653c4762a1bSJed Brown     d_nnz[idx] += 3;
654c4762a1bSJed Brown     d_nnz[idx + 1] += 2;
655c4762a1bSJed Brown     d_nnz[idx + 2] += 2;
656c4762a1bSJed Brown     d_nnz[idx + 3] += 5;
657c4762a1bSJed Brown     d_nnz[idx + 4] += 6;
658c4762a1bSJed Brown     d_nnz[idx + 5] += 6;
659c4762a1bSJed Brown 
660c4762a1bSJed Brown     d_nnz[user->neqs_gen + 2 * gbus[i]] += 3;
661c4762a1bSJed Brown     d_nnz[user->neqs_gen + 2 * gbus[i] + 1] += 3;
662c4762a1bSJed Brown 
663c4762a1bSJed Brown     d_nnz[idx + 6] += 2;
664c4762a1bSJed Brown     d_nnz[idx + 7] += 2;
665c4762a1bSJed Brown     d_nnz[idx + 8] += 5;
666c4762a1bSJed Brown 
667c4762a1bSJed Brown     idx = idx + 9;
668c4762a1bSJed Brown   }
669c4762a1bSJed Brown   start = user->neqs_gen;
670c4762a1bSJed Brown 
671c4762a1bSJed Brown   for (i = 0; i < nbus; i++) {
6729566063dSJacob Faibussowitsch     PetscCall(MatGetRow(user->Ybus, 2 * i, &ncols, NULL, NULL));
673c4762a1bSJed Brown     d_nnz[start + 2 * i] += ncols;
674c4762a1bSJed Brown     d_nnz[start + 2 * i + 1] += ncols;
6759566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(user->Ybus, 2 * i, &ncols, NULL, NULL));
676c4762a1bSJed Brown   }
6779566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(J, 0, d_nnz));
6789566063dSJacob Faibussowitsch   PetscCall(PetscFree(d_nnz));
6793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
680c4762a1bSJed Brown }
681c4762a1bSJed Brown 
682c4762a1bSJed Brown /*
683c4762a1bSJed Brown    J = [df_dx, df_dy
684c4762a1bSJed Brown         dg_dx, dg_dy]
685c4762a1bSJed Brown */
ResidualJacobian(Vec X,Mat J,Mat B,PetscCtx ctx)686*2a8381b2SBarry Smith PetscErrorCode ResidualJacobian(Vec X, Mat J, Mat B, PetscCtx ctx)
687d71ae5a4SJacob Faibussowitsch {
688c4762a1bSJed Brown   Userctx           *user = (Userctx *)ctx;
689c4762a1bSJed Brown   Vec                Xgen, Xnet;
69005e808d2SBarry Smith   const PetscScalar *xgen, *xnet;
691c4762a1bSJed Brown   PetscInt           i, idx = 0;
692c4762a1bSJed Brown   PetscScalar        Vr, Vi, Vm, Vm2;
693c4762a1bSJed Brown   PetscScalar        Eqp, Edp, delta; /* Generator variables */
694c4762a1bSJed Brown   PetscScalar        Efd;
695c4762a1bSJed Brown   PetscScalar        Id, Iq; /* Generator dq axis currents */
696c4762a1bSJed Brown   PetscScalar        Vd, Vq;
697c4762a1bSJed Brown   PetscScalar        val[10];
698c4762a1bSJed Brown   PetscInt           row[2], col[10];
699c4762a1bSJed Brown   PetscInt           net_start = user->neqs_gen;
700c4762a1bSJed Brown   PetscScalar        Zdq_inv[4], det;
701c4762a1bSJed Brown   PetscScalar        dVd_dVr, dVd_dVi, dVq_dVr, dVq_dVi, dVd_ddelta, dVq_ddelta;
702c4762a1bSJed Brown   PetscScalar        dIGr_ddelta, dIGi_ddelta, dIGr_dId, dIGr_dIq, dIGi_dId, dIGi_dIq;
703c4762a1bSJed Brown   PetscScalar        dSE_dEfd;
704c4762a1bSJed Brown   PetscScalar        dVm_dVd, dVm_dVq, dVm_dVr, dVm_dVi;
705c4762a1bSJed Brown   PetscInt           ncols;
706c4762a1bSJed Brown   const PetscInt    *cols;
707c4762a1bSJed Brown   const PetscScalar *yvals;
708c4762a1bSJed Brown   PetscInt           k;
70905e808d2SBarry Smith   PetscScalar        PD, QD, Vm0, Vm4;
71005e808d2SBarry Smith   const PetscScalar *v0;
711c4762a1bSJed Brown   PetscScalar        dPD_dVr, dPD_dVi, dQD_dVr, dQD_dVi;
712c4762a1bSJed Brown   PetscScalar        dIDr_dVr, dIDr_dVi, dIDi_dVr, dIDi_dVi;
713c4762a1bSJed Brown 
714c4762a1bSJed Brown   PetscFunctionBegin;
7159566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(B));
7169566063dSJacob Faibussowitsch   PetscCall(DMCompositeGetLocalVectors(user->dmpgrid, &Xgen, &Xnet));
7179566063dSJacob Faibussowitsch   PetscCall(DMCompositeScatter(user->dmpgrid, X, Xgen, Xnet));
718c4762a1bSJed Brown 
7199566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xgen, &xgen));
7209566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xnet, &xnet));
721c4762a1bSJed Brown 
722c4762a1bSJed Brown   /* Generator subsystem */
723c4762a1bSJed Brown   for (i = 0; i < ngen; i++) {
724c4762a1bSJed Brown     Eqp   = xgen[idx];
725c4762a1bSJed Brown     Edp   = xgen[idx + 1];
726c4762a1bSJed Brown     delta = xgen[idx + 2];
727c4762a1bSJed Brown     Id    = xgen[idx + 4];
728c4762a1bSJed Brown     Iq    = xgen[idx + 5];
729c4762a1bSJed Brown     Efd   = xgen[idx + 6];
730c4762a1bSJed Brown 
731c4762a1bSJed Brown     /*    fgen[idx]   = (-Eqp - (Xd[i] - Xdp[i])*Id + Efd)/Td0p[i]; */
732c4762a1bSJed Brown     row[0] = idx;
7339371c9d4SSatish Balay     col[0] = idx;
7349371c9d4SSatish Balay     col[1] = idx + 4;
7359371c9d4SSatish Balay     col[2] = idx + 6;
7369371c9d4SSatish Balay     val[0] = -1 / Td0p[i];
7379371c9d4SSatish Balay     val[1] = -(Xd[i] - Xdp[i]) / Td0p[i];
7389371c9d4SSatish Balay     val[2] = 1 / Td0p[i];
739c4762a1bSJed Brown 
7409566063dSJacob Faibussowitsch     PetscCall(MatSetValues(J, 1, row, 3, col, val, INSERT_VALUES));
741c4762a1bSJed Brown 
742c4762a1bSJed Brown     /*    fgen[idx+1] = (-Edp + (Xq[i] - Xqp[i])*Iq)/Tq0p[i]; */
743c4762a1bSJed Brown     row[0] = idx + 1;
7449371c9d4SSatish Balay     col[0] = idx + 1;
7459371c9d4SSatish Balay     col[1] = idx + 5;
7469371c9d4SSatish Balay     val[0] = -1 / Tq0p[i];
7479371c9d4SSatish Balay     val[1] = (Xq[i] - Xqp[i]) / Tq0p[i];
7489566063dSJacob Faibussowitsch     PetscCall(MatSetValues(J, 1, row, 2, col, val, INSERT_VALUES));
749c4762a1bSJed Brown 
750c4762a1bSJed Brown     /*    fgen[idx+2] = w - w_s; */
751c4762a1bSJed Brown     row[0] = idx + 2;
7529371c9d4SSatish Balay     col[0] = idx + 2;
7539371c9d4SSatish Balay     col[1] = idx + 3;
7549371c9d4SSatish Balay     val[0] = 0;
7559371c9d4SSatish Balay     val[1] = 1;
7569566063dSJacob Faibussowitsch     PetscCall(MatSetValues(J, 1, row, 2, col, val, INSERT_VALUES));
757c4762a1bSJed Brown 
758c4762a1bSJed Brown     /*    fgen[idx+3] = (TM[i] - Edp*Id - Eqp*Iq - (Xqp[i] - Xdp[i])*Id*Iq - D[i]*(w - w_s))/M[i]; */
759c4762a1bSJed Brown     row[0] = idx + 3;
7609371c9d4SSatish Balay     col[0] = idx;
7619371c9d4SSatish Balay     col[1] = idx + 1;
7629371c9d4SSatish Balay     col[2] = idx + 3;
7639371c9d4SSatish Balay     col[3] = idx + 4;
7649371c9d4SSatish Balay     col[4] = idx + 5;
7659371c9d4SSatish Balay     val[0] = -Iq / M[i];
7669371c9d4SSatish Balay     val[1] = -Id / M[i];
7679371c9d4SSatish Balay     val[2] = -D[i] / M[i];
7689371c9d4SSatish Balay     val[3] = (-Edp - (Xqp[i] - Xdp[i]) * Iq) / M[i];
7699371c9d4SSatish Balay     val[4] = (-Eqp - (Xqp[i] - Xdp[i]) * Id) / M[i];
7709566063dSJacob Faibussowitsch     PetscCall(MatSetValues(J, 1, row, 5, col, val, INSERT_VALUES));
771c4762a1bSJed Brown 
772c4762a1bSJed Brown     Vr = xnet[2 * gbus[i]];     /* Real part of generator terminal voltage */
773c4762a1bSJed Brown     Vi = xnet[2 * gbus[i] + 1]; /* Imaginary part of the generator terminal voltage */
7749566063dSJacob Faibussowitsch     PetscCall(ri2dq(Vr, Vi, delta, &Vd, &Vq));
775c4762a1bSJed Brown 
776c4762a1bSJed Brown     det = Rs[i] * Rs[i] + Xdp[i] * Xqp[i];
777c4762a1bSJed Brown 
778c4762a1bSJed Brown     Zdq_inv[0] = Rs[i] / det;
779c4762a1bSJed Brown     Zdq_inv[1] = Xqp[i] / det;
780c4762a1bSJed Brown     Zdq_inv[2] = -Xdp[i] / det;
781c4762a1bSJed Brown     Zdq_inv[3] = Rs[i] / det;
782c4762a1bSJed Brown 
7839371c9d4SSatish Balay     dVd_dVr    = PetscSinScalar(delta);
7849371c9d4SSatish Balay     dVd_dVi    = -PetscCosScalar(delta);
7859371c9d4SSatish Balay     dVq_dVr    = PetscCosScalar(delta);
7869371c9d4SSatish Balay     dVq_dVi    = PetscSinScalar(delta);
787c4762a1bSJed Brown     dVd_ddelta = Vr * PetscCosScalar(delta) + Vi * PetscSinScalar(delta);
788c4762a1bSJed Brown     dVq_ddelta = -Vr * PetscSinScalar(delta) + Vi * PetscCosScalar(delta);
789c4762a1bSJed Brown 
790c4762a1bSJed Brown     /*    fgen[idx+4] = Zdq_inv[0]*(-Edp + Vd) + Zdq_inv[1]*(-Eqp + Vq) + Id; */
791c4762a1bSJed Brown     row[0] = idx + 4;
7929371c9d4SSatish Balay     col[0] = idx;
7939371c9d4SSatish Balay     col[1] = idx + 1;
7949371c9d4SSatish Balay     col[2] = idx + 2;
7959371c9d4SSatish Balay     val[0] = -Zdq_inv[1];
7969371c9d4SSatish Balay     val[1] = -Zdq_inv[0];
7979371c9d4SSatish Balay     val[2] = Zdq_inv[0] * dVd_ddelta + Zdq_inv[1] * dVq_ddelta;
7989371c9d4SSatish Balay     col[3] = idx + 4;
7999371c9d4SSatish Balay     col[4] = net_start + 2 * gbus[i];
8009371c9d4SSatish Balay     col[5] = net_start + 2 * gbus[i] + 1;
8019371c9d4SSatish Balay     val[3] = 1;
8029371c9d4SSatish Balay     val[4] = Zdq_inv[0] * dVd_dVr + Zdq_inv[1] * dVq_dVr;
8039371c9d4SSatish Balay     val[5] = Zdq_inv[0] * dVd_dVi + Zdq_inv[1] * dVq_dVi;
8049566063dSJacob Faibussowitsch     PetscCall(MatSetValues(J, 1, row, 6, col, val, INSERT_VALUES));
805c4762a1bSJed Brown 
806c4762a1bSJed Brown     /*  fgen[idx+5] = Zdq_inv[2]*(-Edp + Vd) + Zdq_inv[3]*(-Eqp + Vq) + Iq; */
807c4762a1bSJed Brown     row[0] = idx + 5;
8089371c9d4SSatish Balay     col[0] = idx;
8099371c9d4SSatish Balay     col[1] = idx + 1;
8109371c9d4SSatish Balay     col[2] = idx + 2;
8119371c9d4SSatish Balay     val[0] = -Zdq_inv[3];
8129371c9d4SSatish Balay     val[1] = -Zdq_inv[2];
8139371c9d4SSatish Balay     val[2] = Zdq_inv[2] * dVd_ddelta + Zdq_inv[3] * dVq_ddelta;
8149371c9d4SSatish Balay     col[3] = idx + 5;
8159371c9d4SSatish Balay     col[4] = net_start + 2 * gbus[i];
8169371c9d4SSatish Balay     col[5] = net_start + 2 * gbus[i] + 1;
8179371c9d4SSatish Balay     val[3] = 1;
8189371c9d4SSatish Balay     val[4] = Zdq_inv[2] * dVd_dVr + Zdq_inv[3] * dVq_dVr;
8199371c9d4SSatish Balay     val[5] = Zdq_inv[2] * dVd_dVi + Zdq_inv[3] * dVq_dVi;
8209566063dSJacob Faibussowitsch     PetscCall(MatSetValues(J, 1, row, 6, col, val, INSERT_VALUES));
821c4762a1bSJed Brown 
822c4762a1bSJed Brown     dIGr_ddelta = Id * PetscCosScalar(delta) - Iq * PetscSinScalar(delta);
823c4762a1bSJed Brown     dIGi_ddelta = Id * PetscSinScalar(delta) + Iq * PetscCosScalar(delta);
8249371c9d4SSatish Balay     dIGr_dId    = PetscSinScalar(delta);
8259371c9d4SSatish Balay     dIGr_dIq    = PetscCosScalar(delta);
8269371c9d4SSatish Balay     dIGi_dId    = -PetscCosScalar(delta);
8279371c9d4SSatish Balay     dIGi_dIq    = PetscSinScalar(delta);
828c4762a1bSJed Brown 
829c4762a1bSJed Brown     /* fnet[2*gbus[i]]   -= IGi; */
830c4762a1bSJed Brown     row[0] = net_start + 2 * gbus[i];
8319371c9d4SSatish Balay     col[0] = idx + 2;
8329371c9d4SSatish Balay     col[1] = idx + 4;
8339371c9d4SSatish Balay     col[2] = idx + 5;
8349371c9d4SSatish Balay     val[0] = -dIGi_ddelta;
8359371c9d4SSatish Balay     val[1] = -dIGi_dId;
8369371c9d4SSatish Balay     val[2] = -dIGi_dIq;
8379566063dSJacob Faibussowitsch     PetscCall(MatSetValues(J, 1, row, 3, col, val, INSERT_VALUES));
838c4762a1bSJed Brown 
839c4762a1bSJed Brown     /* fnet[2*gbus[i]+1]   -= IGr; */
840c4762a1bSJed Brown     row[0] = net_start + 2 * gbus[i] + 1;
8419371c9d4SSatish Balay     col[0] = idx + 2;
8429371c9d4SSatish Balay     col[1] = idx + 4;
8439371c9d4SSatish Balay     col[2] = idx + 5;
8449371c9d4SSatish Balay     val[0] = -dIGr_ddelta;
8459371c9d4SSatish Balay     val[1] = -dIGr_dId;
8469371c9d4SSatish Balay     val[2] = -dIGr_dIq;
8479566063dSJacob Faibussowitsch     PetscCall(MatSetValues(J, 1, row, 3, col, val, INSERT_VALUES));
848c4762a1bSJed Brown 
849c4762a1bSJed Brown     Vm = PetscSqrtScalar(Vd * Vd + Vq * Vq);
850c4762a1bSJed Brown 
851c4762a1bSJed Brown     /*    fgen[idx+6] = (-KE[i]*Efd - SE + VR)/TE[i]; */
852c4762a1bSJed Brown     /*    SE  = k1[i]*PetscExpScalar(k2[i]*Efd); */
853c4762a1bSJed Brown 
854c4762a1bSJed Brown     dSE_dEfd = k1[i] * k2[i] * PetscExpScalar(k2[i] * Efd);
855c4762a1bSJed Brown 
856c4762a1bSJed Brown     row[0] = idx + 6;
8579371c9d4SSatish Balay     col[0] = idx + 6;
8589371c9d4SSatish Balay     col[1] = idx + 8;
8599371c9d4SSatish Balay     val[0] = (-KE[i] - dSE_dEfd) / TE[i];
8609371c9d4SSatish Balay     val[1] = 1 / TE[i];
8619566063dSJacob Faibussowitsch     PetscCall(MatSetValues(J, 1, row, 2, col, val, INSERT_VALUES));
862c4762a1bSJed Brown 
863c4762a1bSJed Brown     /* Exciter differential equations */
864c4762a1bSJed Brown 
865c4762a1bSJed Brown     /*    fgen[idx+7] = (-RF + KF[i]*Efd/TF[i])/TF[i]; */
866c4762a1bSJed Brown     row[0] = idx + 7;
8679371c9d4SSatish Balay     col[0] = idx + 6;
8689371c9d4SSatish Balay     col[1] = idx + 7;
8699371c9d4SSatish Balay     val[0] = (KF[i] / TF[i]) / TF[i];
8709371c9d4SSatish Balay     val[1] = -1 / TF[i];
8719566063dSJacob Faibussowitsch     PetscCall(MatSetValues(J, 1, row, 2, col, val, INSERT_VALUES));
872c4762a1bSJed Brown 
873c4762a1bSJed Brown     /*    fgen[idx+8] = (-VR + KA[i]*RF - KA[i]*KF[i]*Efd/TF[i] + KA[i]*(Vref[i] - Vm))/TA[i]; */
874c4762a1bSJed Brown     /* Vm = (Vd^2 + Vq^2)^0.5; */
875c4762a1bSJed Brown 
876c4762a1bSJed Brown     row[0] = idx + 8;
877c4762a1bSJed Brown     if (VRatmax[i]) {
8789371c9d4SSatish Balay       col[0] = idx + 8;
8799371c9d4SSatish Balay       val[0] = 1.0;
8809566063dSJacob Faibussowitsch       PetscCall(MatSetValues(J, 1, row, 1, col, val, INSERT_VALUES));
881c4762a1bSJed Brown     } else if (VRatmin[i]) {
8829371c9d4SSatish Balay       col[0] = idx + 8;
8839371c9d4SSatish Balay       val[0] = -1.0;
8849566063dSJacob Faibussowitsch       PetscCall(MatSetValues(J, 1, row, 1, col, val, INSERT_VALUES));
885c4762a1bSJed Brown     } else {
8869371c9d4SSatish Balay       dVm_dVd = Vd / Vm;
8879371c9d4SSatish Balay       dVm_dVq = Vq / Vm;
888c4762a1bSJed Brown       dVm_dVr = dVm_dVd * dVd_dVr + dVm_dVq * dVq_dVr;
889c4762a1bSJed Brown       dVm_dVi = dVm_dVd * dVd_dVi + dVm_dVq * dVq_dVi;
890c4762a1bSJed Brown       row[0]  = idx + 8;
8919371c9d4SSatish Balay       col[0]  = idx + 6;
8929371c9d4SSatish Balay       col[1]  = idx + 7;
8939371c9d4SSatish Balay       col[2]  = idx + 8;
8949371c9d4SSatish Balay       val[0]  = -(KA[i] * KF[i] / TF[i]) / TA[i];
8959371c9d4SSatish Balay       val[1]  = KA[i] / TA[i];
8969371c9d4SSatish Balay       val[2]  = -1 / TA[i];
8979371c9d4SSatish Balay       col[3]  = net_start + 2 * gbus[i];
8989371c9d4SSatish Balay       col[4]  = net_start + 2 * gbus[i] + 1;
8999371c9d4SSatish Balay       val[3]  = -KA[i] * dVm_dVr / TA[i];
9009371c9d4SSatish Balay       val[4]  = -KA[i] * dVm_dVi / TA[i];
9019566063dSJacob Faibussowitsch       PetscCall(MatSetValues(J, 1, row, 5, col, val, INSERT_VALUES));
902c4762a1bSJed Brown     }
903c4762a1bSJed Brown     idx = idx + 9;
904c4762a1bSJed Brown   }
905c4762a1bSJed Brown 
906c4762a1bSJed Brown   for (i = 0; i < nbus; i++) {
9079566063dSJacob Faibussowitsch     PetscCall(MatGetRow(user->Ybus, 2 * i, &ncols, &cols, &yvals));
908c4762a1bSJed Brown     row[0] = net_start + 2 * i;
909c4762a1bSJed Brown     for (k = 0; k < ncols; k++) {
910c4762a1bSJed Brown       col[k] = net_start + cols[k];
911c4762a1bSJed Brown       val[k] = yvals[k];
912c4762a1bSJed Brown     }
9139566063dSJacob Faibussowitsch     PetscCall(MatSetValues(J, 1, row, ncols, col, val, INSERT_VALUES));
9149566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(user->Ybus, 2 * i, &ncols, &cols, &yvals));
915c4762a1bSJed Brown 
9169566063dSJacob Faibussowitsch     PetscCall(MatGetRow(user->Ybus, 2 * i + 1, &ncols, &cols, &yvals));
917c4762a1bSJed Brown     row[0] = net_start + 2 * i + 1;
918c4762a1bSJed Brown     for (k = 0; k < ncols; k++) {
919c4762a1bSJed Brown       col[k] = net_start + cols[k];
920c4762a1bSJed Brown       val[k] = yvals[k];
921c4762a1bSJed Brown     }
9229566063dSJacob Faibussowitsch     PetscCall(MatSetValues(J, 1, row, ncols, col, val, INSERT_VALUES));
9239566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(user->Ybus, 2 * i + 1, &ncols, &cols, &yvals));
924c4762a1bSJed Brown   }
925c4762a1bSJed Brown 
9269566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(J, MAT_FLUSH_ASSEMBLY));
9279566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(J, MAT_FLUSH_ASSEMBLY));
928c4762a1bSJed Brown 
9299566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(user->V0, &v0));
930c4762a1bSJed Brown   for (i = 0; i < nload; i++) {
931c4762a1bSJed Brown     Vr  = xnet[2 * lbus[i]];     /* Real part of load bus voltage */
932c4762a1bSJed Brown     Vi  = xnet[2 * lbus[i] + 1]; /* Imaginary part of the load bus voltage */
9339371c9d4SSatish Balay     Vm  = PetscSqrtScalar(Vr * Vr + Vi * Vi);
9349371c9d4SSatish Balay     Vm2 = Vm * Vm;
9359371c9d4SSatish Balay     Vm4 = Vm2 * Vm2;
936c4762a1bSJed Brown     Vm0 = PetscSqrtScalar(v0[2 * lbus[i]] * v0[2 * lbus[i]] + v0[2 * lbus[i] + 1] * v0[2 * lbus[i] + 1]);
937c4762a1bSJed Brown     PD = QD = 0.0;
938c4762a1bSJed Brown     dPD_dVr = dPD_dVi = dQD_dVr = dQD_dVi = 0.0;
939c4762a1bSJed Brown     for (k = 0; k < ld_nsegsp[i]; k++) {
94057508eceSPierre Jolivet       PD += ld_alphap[k] * PD0[i] * PetscPowScalar(Vm / Vm0, ld_betap[k]);
94157508eceSPierre Jolivet       dPD_dVr += ld_alphap[k] * ld_betap[k] * PD0[i] * PetscPowScalar(1 / Vm0, ld_betap[k]) * Vr * PetscPowScalar(Vm, ld_betap[k] - 2);
94257508eceSPierre Jolivet       dPD_dVi += ld_alphap[k] * ld_betap[k] * PD0[i] * PetscPowScalar(1 / Vm0, ld_betap[k]) * Vi * PetscPowScalar(Vm, ld_betap[k] - 2);
943c4762a1bSJed Brown     }
944c4762a1bSJed Brown     for (k = 0; k < ld_nsegsq[i]; k++) {
94557508eceSPierre Jolivet       QD += ld_alphaq[k] * QD0[i] * PetscPowScalar(Vm / Vm0, ld_betaq[k]);
94657508eceSPierre Jolivet       dQD_dVr += ld_alphaq[k] * ld_betaq[k] * QD0[i] * PetscPowScalar(1 / Vm0, ld_betaq[k]) * Vr * PetscPowScalar(Vm, ld_betaq[k] - 2);
94757508eceSPierre Jolivet       dQD_dVi += ld_alphaq[k] * ld_betaq[k] * QD0[i] * PetscPowScalar(1 / Vm0, ld_betaq[k]) * Vi * PetscPowScalar(Vm, ld_betaq[k] - 2);
948c4762a1bSJed Brown     }
949c4762a1bSJed Brown 
950c4762a1bSJed Brown     /*    IDr = (PD*Vr + QD*Vi)/Vm2; */
951c4762a1bSJed Brown     /*    IDi = (-QD*Vr + PD*Vi)/Vm2; */
952c4762a1bSJed Brown 
953c4762a1bSJed Brown     dIDr_dVr = (dPD_dVr * Vr + dQD_dVr * Vi + PD) / Vm2 - ((PD * Vr + QD * Vi) * 2 * Vr) / Vm4;
954c4762a1bSJed Brown     dIDr_dVi = (dPD_dVi * Vr + dQD_dVi * Vi + QD) / Vm2 - ((PD * Vr + QD * Vi) * 2 * Vi) / Vm4;
955c4762a1bSJed Brown 
956c4762a1bSJed Brown     dIDi_dVr = (-dQD_dVr * Vr + dPD_dVr * Vi - QD) / Vm2 - ((-QD * Vr + PD * Vi) * 2 * Vr) / Vm4;
957c4762a1bSJed Brown     dIDi_dVi = (-dQD_dVi * Vr + dPD_dVi * Vi + PD) / Vm2 - ((-QD * Vr + PD * Vi) * 2 * Vi) / Vm4;
958c4762a1bSJed Brown 
959c4762a1bSJed Brown     /*    fnet[2*lbus[i]]   += IDi; */
960c4762a1bSJed Brown     row[0] = net_start + 2 * lbus[i];
9619371c9d4SSatish Balay     col[0] = net_start + 2 * lbus[i];
9629371c9d4SSatish Balay     col[1] = net_start + 2 * lbus[i] + 1;
9639371c9d4SSatish Balay     val[0] = dIDi_dVr;
9649371c9d4SSatish Balay     val[1] = dIDi_dVi;
9659566063dSJacob Faibussowitsch     PetscCall(MatSetValues(J, 1, row, 2, col, val, ADD_VALUES));
966c4762a1bSJed Brown     /*    fnet[2*lbus[i]+1] += IDr; */
967c4762a1bSJed Brown     row[0] = net_start + 2 * lbus[i] + 1;
9689371c9d4SSatish Balay     col[0] = net_start + 2 * lbus[i];
9699371c9d4SSatish Balay     col[1] = net_start + 2 * lbus[i] + 1;
9709371c9d4SSatish Balay     val[0] = dIDr_dVr;
9719371c9d4SSatish Balay     val[1] = dIDr_dVi;
9729566063dSJacob Faibussowitsch     PetscCall(MatSetValues(J, 1, row, 2, col, val, ADD_VALUES));
973c4762a1bSJed Brown   }
9749566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(user->V0, &v0));
975c4762a1bSJed Brown 
9769566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xgen, &xgen));
9779566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xnet, &xnet));
978c4762a1bSJed Brown 
9799566063dSJacob Faibussowitsch   PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid, &Xgen, &Xnet));
980c4762a1bSJed Brown 
9819566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
9829566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
9833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
984c4762a1bSJed Brown }
985c4762a1bSJed Brown 
986c4762a1bSJed Brown /*
987c4762a1bSJed Brown    J = [I, 0
988c4762a1bSJed Brown         dg_dx, dg_dy]
989c4762a1bSJed Brown */
AlgJacobian(SNES snes,Vec X,Mat A,Mat B,PetscCtx ctx)990*2a8381b2SBarry Smith PetscErrorCode AlgJacobian(SNES snes, Vec X, Mat A, Mat B, PetscCtx ctx)
991d71ae5a4SJacob Faibussowitsch {
992c4762a1bSJed Brown   Userctx *user = (Userctx *)ctx;
993c4762a1bSJed Brown 
994c4762a1bSJed Brown   PetscFunctionBegin;
9959566063dSJacob Faibussowitsch   PetscCall(ResidualJacobian(X, A, B, ctx));
9969566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A, MAT_KEEP_NONZERO_PATTERN, PETSC_TRUE));
9979566063dSJacob Faibussowitsch   PetscCall(MatZeroRowsIS(A, user->is_diff, 1.0, NULL, NULL));
9983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
999c4762a1bSJed Brown }
1000c4762a1bSJed Brown 
1001c4762a1bSJed Brown /*
1002c4762a1bSJed Brown    J = [-df_dx, -df_dy
1003c4762a1bSJed Brown         dg_dx, dg_dy]
1004c4762a1bSJed Brown */
1005c4762a1bSJed Brown 
RHSJacobian(TS ts,PetscReal t,Vec X,Mat A,Mat B,PetscCtx ctx)1006*2a8381b2SBarry Smith PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec X, Mat A, Mat B, PetscCtx ctx)
1007d71ae5a4SJacob Faibussowitsch {
1008c4762a1bSJed Brown   Userctx *user = (Userctx *)ctx;
1009c4762a1bSJed Brown 
1010c4762a1bSJed Brown   PetscFunctionBegin;
1011c4762a1bSJed Brown   user->t = t;
1012c4762a1bSJed Brown 
10139566063dSJacob Faibussowitsch   PetscCall(ResidualJacobian(X, A, B, user));
10143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1015c4762a1bSJed Brown }
1016c4762a1bSJed Brown 
1017c4762a1bSJed Brown /*
1018c4762a1bSJed Brown    J = [df_dx-aI, df_dy
1019c4762a1bSJed Brown         dg_dx, dg_dy]
1020c4762a1bSJed Brown */
1021c4762a1bSJed Brown 
IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,Userctx * user)1022d71ae5a4SJacob Faibussowitsch PetscErrorCode IJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat A, Mat B, Userctx *user)
1023d71ae5a4SJacob Faibussowitsch {
1024c4762a1bSJed Brown   PetscScalar atmp = (PetscScalar)a;
1025c4762a1bSJed Brown   PetscInt    i, row;
1026c4762a1bSJed Brown 
1027c4762a1bSJed Brown   PetscFunctionBegin;
1028c4762a1bSJed Brown   user->t = t;
1029c4762a1bSJed Brown 
10309566063dSJacob Faibussowitsch   PetscCall(RHSJacobian(ts, t, X, A, B, user));
10319566063dSJacob Faibussowitsch   PetscCall(MatScale(B, -1.0));
1032c4762a1bSJed Brown   for (i = 0; i < ngen; i++) {
1033c4762a1bSJed Brown     row = 9 * i;
10349566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES));
1035c4762a1bSJed Brown     row = 9 * i + 1;
10369566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES));
1037c4762a1bSJed Brown     row = 9 * i + 2;
10389566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES));
1039c4762a1bSJed Brown     row = 9 * i + 3;
10409566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES));
1041c4762a1bSJed Brown     row = 9 * i + 6;
10429566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES));
1043c4762a1bSJed Brown     row = 9 * i + 7;
10449566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES));
1045c4762a1bSJed Brown     row = 9 * i + 8;
10469566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES));
1047c4762a1bSJed Brown   }
10489566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
10499566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
10503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1051c4762a1bSJed Brown }
1052c4762a1bSJed Brown 
main(int argc,char ** argv)1053d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
1054d71ae5a4SJacob Faibussowitsch {
1055c4762a1bSJed Brown   TS                 ts;
1056c4762a1bSJed Brown   SNES               snes_alg;
1057c4762a1bSJed Brown   PetscMPIInt        size;
1058c4762a1bSJed Brown   Userctx            user;
1059c4762a1bSJed Brown   PetscViewer        Xview, Ybusview, viewer;
1060c4762a1bSJed Brown   Vec                X, F_alg;
1061c4762a1bSJed Brown   Mat                J, A;
1062c4762a1bSJed Brown   PetscInt           i, idx, *idx2;
1063c4762a1bSJed Brown   Vec                Xdot;
1064c4762a1bSJed Brown   PetscScalar       *x, *mat, *amat;
1065c4762a1bSJed Brown   const PetscScalar *rmat;
1066c4762a1bSJed Brown   Vec                vatol;
1067c4762a1bSJed Brown   PetscInt          *direction;
1068c4762a1bSJed Brown   PetscBool         *terminate;
1069c4762a1bSJed Brown   const PetscInt    *idx3;
1070c4762a1bSJed Brown   PetscScalar       *vatoli;
1071c4762a1bSJed Brown   PetscInt           k;
1072c4762a1bSJed Brown 
1073327415f7SBarry Smith   PetscFunctionBeginUser;
10749566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, "petscoptions", help));
10759566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
10763c633725SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs");
1077c4762a1bSJed Brown 
1078c4762a1bSJed Brown   user.neqs_gen   = 9 * ngen; /* # eqs. for generator subsystem */
1079c4762a1bSJed Brown   user.neqs_net   = 2 * nbus; /* # eqs. for network subsystem   */
1080c4762a1bSJed Brown   user.neqs_pgrid = user.neqs_gen + user.neqs_net;
1081c4762a1bSJed Brown 
1082c4762a1bSJed Brown   /* Create indices for differential and algebraic equations */
1083c4762a1bSJed Brown 
10849566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(7 * ngen, &idx2));
1085c4762a1bSJed Brown   for (i = 0; i < ngen; i++) {
10869371c9d4SSatish Balay     idx2[7 * i]     = 9 * i;
10879371c9d4SSatish Balay     idx2[7 * i + 1] = 9 * i + 1;
10889371c9d4SSatish Balay     idx2[7 * i + 2] = 9 * i + 2;
10899371c9d4SSatish Balay     idx2[7 * i + 3] = 9 * i + 3;
10909371c9d4SSatish Balay     idx2[7 * i + 4] = 9 * i + 6;
10919371c9d4SSatish Balay     idx2[7 * i + 5] = 9 * i + 7;
10929371c9d4SSatish Balay     idx2[7 * i + 6] = 9 * i + 8;
1093c4762a1bSJed Brown   }
10949566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, 7 * ngen, idx2, PETSC_COPY_VALUES, &user.is_diff));
10959566063dSJacob Faibussowitsch   PetscCall(ISComplement(user.is_diff, 0, user.neqs_pgrid, &user.is_alg));
10969566063dSJacob Faibussowitsch   PetscCall(PetscFree(idx2));
1097c4762a1bSJed Brown 
1098c4762a1bSJed Brown   /* Read initial voltage vector and Ybus */
10999566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "X.bin", FILE_MODE_READ, &Xview));
11009566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "Ybus.bin", FILE_MODE_READ, &Ybusview));
1101c4762a1bSJed Brown 
11029566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &user.V0));
11039566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(user.V0, PETSC_DECIDE, user.neqs_net));
11049566063dSJacob Faibussowitsch   PetscCall(VecLoad(user.V0, Xview));
1105c4762a1bSJed Brown 
11069566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &user.Ybus));
11079566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(user.Ybus, PETSC_DECIDE, PETSC_DECIDE, user.neqs_net, user.neqs_net));
11089566063dSJacob Faibussowitsch   PetscCall(MatSetType(user.Ybus, MATBAIJ));
11099566063dSJacob Faibussowitsch   /*  PetscCall(MatSetBlockSize(user.Ybus,2)); */
11109566063dSJacob Faibussowitsch   PetscCall(MatLoad(user.Ybus, Ybusview));
1111c4762a1bSJed Brown 
1112c4762a1bSJed Brown   /* Set run time options */
1113d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Transient stability fault options", "");
1114c4762a1bSJed Brown   {
1115c4762a1bSJed Brown     user.tfaulton     = 1.0;
1116c4762a1bSJed Brown     user.tfaultoff    = 1.2;
1117c4762a1bSJed Brown     user.Rfault       = 0.0001;
1118c4762a1bSJed Brown     user.setisdiff    = PETSC_FALSE;
1119c4762a1bSJed Brown     user.semiexplicit = PETSC_FALSE;
1120c4762a1bSJed Brown     user.faultbus     = 8;
11219566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-tfaulton", "", "", user.tfaulton, &user.tfaulton, NULL));
11229566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-tfaultoff", "", "", user.tfaultoff, &user.tfaultoff, NULL));
11239566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-faultbus", "", "", user.faultbus, &user.faultbus, NULL));
1124c4762a1bSJed Brown     user.t0   = 0.0;
1125c4762a1bSJed Brown     user.tmax = 5.0;
11269566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-t0", "", "", user.t0, &user.t0, NULL));
11279566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-tmax", "", "", user.tmax, &user.tmax, NULL));
11289566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-setisdiff", "", "", user.setisdiff, &user.setisdiff, NULL));
11299566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-dae_semiexplicit", "", "", user.semiexplicit, &user.semiexplicit, NULL));
1130c4762a1bSJed Brown   }
1131d0609cedSBarry Smith   PetscOptionsEnd();
1132c4762a1bSJed Brown 
11339566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&Xview));
11349566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&Ybusview));
1135c4762a1bSJed Brown 
1136c4762a1bSJed Brown   /* Create DMs for generator and network subsystems */
11379566063dSJacob Faibussowitsch   PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, user.neqs_gen, 1, 1, NULL, &user.dmgen));
11389566063dSJacob Faibussowitsch   PetscCall(DMSetOptionsPrefix(user.dmgen, "dmgen_"));
11399566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(user.dmgen));
11409566063dSJacob Faibussowitsch   PetscCall(DMSetUp(user.dmgen));
11419566063dSJacob Faibussowitsch   PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, user.neqs_net, 1, 1, NULL, &user.dmnet));
11429566063dSJacob Faibussowitsch   PetscCall(DMSetOptionsPrefix(user.dmnet, "dmnet_"));
11439566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(user.dmnet));
11449566063dSJacob Faibussowitsch   PetscCall(DMSetUp(user.dmnet));
1145c4762a1bSJed Brown   /* Create a composite DM packer and add the two DMs */
11469566063dSJacob Faibussowitsch   PetscCall(DMCompositeCreate(PETSC_COMM_WORLD, &user.dmpgrid));
11479566063dSJacob Faibussowitsch   PetscCall(DMSetOptionsPrefix(user.dmpgrid, "pgrid_"));
11489566063dSJacob Faibussowitsch   PetscCall(DMCompositeAddDM(user.dmpgrid, user.dmgen));
11499566063dSJacob Faibussowitsch   PetscCall(DMCompositeAddDM(user.dmpgrid, user.dmnet));
1150c4762a1bSJed Brown 
11519566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(user.dmpgrid, &X));
1152c4762a1bSJed Brown 
11539566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &J));
11549566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(J, PETSC_DECIDE, PETSC_DECIDE, user.neqs_pgrid, user.neqs_pgrid));
11559566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(J));
11569566063dSJacob Faibussowitsch   PetscCall(PreallocateJacobian(J, &user));
1157c4762a1bSJed Brown 
1158c4762a1bSJed Brown   /* Create matrix to save solutions at each time step */
1159c4762a1bSJed Brown   user.stepnum = 0;
1160c4762a1bSJed Brown 
11619566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, user.neqs_pgrid + 1, 1002, NULL, &user.Sol));
1162c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1163c4762a1bSJed Brown      Create timestepping solver context
1164c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
11659566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
11669566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
1167c4762a1bSJed Brown   if (user.semiexplicit) {
11689566063dSJacob Faibussowitsch     PetscCall(TSSetType(ts, TSRK));
11699566063dSJacob Faibussowitsch     PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, &user));
11709566063dSJacob Faibussowitsch     PetscCall(TSSetRHSJacobian(ts, J, J, RHSJacobian, &user));
1171c4762a1bSJed Brown   } else {
11729566063dSJacob Faibussowitsch     PetscCall(TSSetType(ts, TSCN));
11739566063dSJacob Faibussowitsch     PetscCall(TSSetEquationType(ts, TS_EQ_DAE_IMPLICIT_INDEX1));
11748434afd1SBarry Smith     PetscCall(TSSetIFunction(ts, NULL, (TSIFunctionFn *)IFunction, &user));
11758434afd1SBarry Smith     PetscCall(TSSetIJacobian(ts, J, J, (TSIJacobianFn *)IJacobian, &user));
1176c4762a1bSJed Brown   }
11779566063dSJacob Faibussowitsch   PetscCall(TSSetApplicationContext(ts, &user));
1178c4762a1bSJed Brown 
1179c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1180c4762a1bSJed Brown      Set initial conditions
1181c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
11829566063dSJacob Faibussowitsch   PetscCall(SetInitialGuess(X, &user));
1183c4762a1bSJed Brown   /* Just to set up the Jacobian structure */
1184c4762a1bSJed Brown 
11859566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(X, &Xdot));
11869566063dSJacob Faibussowitsch   PetscCall(IJacobian(ts, 0.0, X, Xdot, 0.0, J, J, &user));
11879566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Xdot));
1188c4762a1bSJed Brown 
1189c4762a1bSJed Brown   /* Save initial solution */
1190c4762a1bSJed Brown 
1191c4762a1bSJed Brown   idx = user.stepnum * (user.neqs_pgrid + 1);
11929566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArray(user.Sol, &mat));
11939566063dSJacob Faibussowitsch   PetscCall(VecGetArray(X, &x));
1194c4762a1bSJed Brown 
1195c4762a1bSJed Brown   mat[idx] = 0.0;
1196c4762a1bSJed Brown 
11979566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(mat + idx + 1, x, user.neqs_pgrid));
11989566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArray(user.Sol, &mat));
11999566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(X, &x));
1200c4762a1bSJed Brown   user.stepnum++;
1201c4762a1bSJed Brown 
12029566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, user.tmax));
12039566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
12049566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, 0.01));
12059566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
12069566063dSJacob Faibussowitsch   PetscCall(TSSetPostStep(ts, SaveSolution));
12079566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, X));
1208c4762a1bSJed Brown 
120957508eceSPierre Jolivet   PetscCall(PetscMalloc1(2 * ngen + 2, &direction));
121057508eceSPierre Jolivet   PetscCall(PetscMalloc1(2 * ngen + 2, &terminate));
1211c4762a1bSJed Brown   direction[0] = direction[1] = 1;
1212c4762a1bSJed Brown   terminate[0] = terminate[1] = PETSC_FALSE;
1213c4762a1bSJed Brown   for (i = 0; i < ngen; i++) {
12149371c9d4SSatish Balay     direction[2 + 2 * i]     = -1;
12159371c9d4SSatish Balay     direction[2 + 2 * i + 1] = 1;
1216c4762a1bSJed Brown     terminate[2 + 2 * i] = terminate[2 + 2 * i + 1] = PETSC_FALSE;
1217c4762a1bSJed Brown   }
1218c4762a1bSJed Brown 
12199566063dSJacob Faibussowitsch   PetscCall(TSSetEventHandler(ts, 2 * ngen + 2, direction, terminate, EventFunction, PostEventFunction, (void *)&user));
1220c4762a1bSJed Brown 
1221c4762a1bSJed Brown   if (user.semiexplicit) {
1222c4762a1bSJed Brown     /* Use a semi-explicit approach with the time-stepping done by an explicit method and the
1223c4762a1bSJed Brown        algrebraic part solved via PostStage and PostEvaluate callbacks
1224c4762a1bSJed Brown     */
12259566063dSJacob Faibussowitsch     PetscCall(TSSetType(ts, TSRK));
12269566063dSJacob Faibussowitsch     PetscCall(TSSetPostStage(ts, PostStage));
12279566063dSJacob Faibussowitsch     PetscCall(TSSetPostEvaluate(ts, PostEvaluate));
1228c4762a1bSJed Brown   }
1229c4762a1bSJed Brown 
1230c4762a1bSJed Brown   if (user.setisdiff) {
1231c4762a1bSJed Brown     /* Create vector of absolute tolerances and set the algebraic part to infinity */
12329566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(X, &vatol));
12339566063dSJacob Faibussowitsch     PetscCall(VecSet(vatol, 100000.0));
12349566063dSJacob Faibussowitsch     PetscCall(VecGetArray(vatol, &vatoli));
12359566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(user.is_diff, &idx3));
1236c4762a1bSJed Brown     for (k = 0; k < 7 * ngen; k++) vatoli[idx3[k]] = 1e-2;
12379566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(vatol, &vatoli));
1238c4762a1bSJed Brown   }
1239c4762a1bSJed Brown 
1240c4762a1bSJed Brown   /* Create the nonlinear solver for solving the algebraic system */
1241c4762a1bSJed Brown   /* Note that although the algebraic system needs to be solved only for
1242c4762a1bSJed Brown      Idq and V, we reuse the entire system including xgen. The xgen
1243c4762a1bSJed Brown      variables are held constant by setting their residuals to 0 and
1244c4762a1bSJed Brown      putting a 1 on the Jacobian diagonal for xgen rows
1245c4762a1bSJed Brown   */
1246c4762a1bSJed Brown 
12479566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(X, &F_alg));
12489566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes_alg));
12499566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes_alg, F_alg, AlgFunction, &user));
12509566063dSJacob Faibussowitsch   PetscCall(SNESSetJacobian(snes_alg, J, J, AlgJacobian, &user));
12519566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes_alg));
1252c4762a1bSJed Brown 
1253c4762a1bSJed Brown   user.snes_alg = snes_alg;
1254c4762a1bSJed Brown 
1255c4762a1bSJed Brown   /* Solve */
12569566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, X));
1257c4762a1bSJed Brown 
12589566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(user.Sol, MAT_FINAL_ASSEMBLY));
12599566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(user.Sol, MAT_FINAL_ASSEMBLY));
1260c4762a1bSJed Brown 
12619566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, user.neqs_pgrid + 1, user.stepnum, NULL, &A));
12629566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArrayRead(user.Sol, &rmat));
12639566063dSJacob Faibussowitsch   PetscCall(MatDenseGetArray(A, &amat));
12649566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(amat, rmat, user.stepnum * (user.neqs_pgrid + 1)));
12659566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArray(A, &amat));
12669566063dSJacob Faibussowitsch   PetscCall(MatDenseRestoreArrayRead(user.Sol, &rmat));
12679566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_SELF, "out.bin", FILE_MODE_WRITE, &viewer));
12689566063dSJacob Faibussowitsch   PetscCall(MatView(A, viewer));
12699566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&viewer));
12709566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
1271c4762a1bSJed Brown 
12729566063dSJacob Faibussowitsch   PetscCall(PetscFree(direction));
12739566063dSJacob Faibussowitsch   PetscCall(PetscFree(terminate));
12749566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes_alg));
12759566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&F_alg));
12769566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
12779566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user.Ybus));
12789566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user.Sol));
12799566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X));
12809566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.V0));
12819566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&user.dmgen));
12829566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&user.dmnet));
12839566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&user.dmpgrid));
12849566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&user.is_diff));
12859566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&user.is_alg));
12869566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
128748a46eb9SPierre Jolivet   if (user.setisdiff) PetscCall(VecDestroy(&vatol));
12889566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
1289b122ec5aSJacob Faibussowitsch   return 0;
1290c4762a1bSJed Brown }
1291c4762a1bSJed Brown 
1292c4762a1bSJed Brown /*TEST
1293c4762a1bSJed Brown 
1294c4762a1bSJed Brown    build:
1295dfd57a17SPierre Jolivet       requires: double !complex !defined(PETSC_USE_64BIT_INDICES)
1296c4762a1bSJed Brown 
1297c4762a1bSJed Brown    test:
1298c4762a1bSJed Brown       suffix: implicit
1299fe4ad979SIlya Fursov       args: -ts_monitor -snes_monitor_short
1300c4762a1bSJed Brown       localrunfiles: petscoptions X.bin Ybus.bin
1301c4762a1bSJed Brown 
1302c4762a1bSJed Brown    test:
1303c4762a1bSJed Brown       suffix: semiexplicit
1304c61711c8SStefano Zampini       args: -ts_monitor -dae_semiexplicit -snes_error_if_not_converged -ts_rk_type 2a
1305c4762a1bSJed Brown       localrunfiles: petscoptions X.bin Ybus.bin
1306c4762a1bSJed Brown 
1307c4762a1bSJed Brown    test:
1308c4762a1bSJed Brown       suffix: steprestart
130995a2cb33SBarry Smith       # needs ARKIMEX methods with all implicit stages since the mass matrix is not the identity
1310fe4ad979SIlya Fursov       args: -ts_monitor -snes_monitor_short -ts_type arkimex -ts_arkimex_type prssp2
1311c4762a1bSJed Brown       localrunfiles: petscoptions X.bin Ybus.bin
1312c4762a1bSJed Brown 
1313c4762a1bSJed Brown TEST*/
1314