| Line | Branch | Exec | Source |
|---|---|---|---|
| 1 | /* | ||
| 2 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 3 | SLEPc - Scalable Library for Eigenvalue Problem Computations | ||
| 4 | Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain | ||
| 5 | |||
| 6 | This file is part of SLEPc. | ||
| 7 | SLEPc is distributed under a 2-clause BSD license (see LICENSE). | ||
| 8 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 9 | */ | ||
| 10 | |||
| 11 | static char help[] = "Power grid small signal stability analysis of WECC 9 bus system.\n\n\ | ||
| 12 | This example is based on the 9-bus (node) example given in the book Power\n\ | ||
| 13 | Systems Dynamics and Stability (Chapter 8) by P. Sauer and M. A. Pai.\n\ | ||
| 14 | The power grid in this example consists of 9 buses (nodes), 3 generators,\n\ | ||
| 15 | 3 loads, and 9 transmission lines. The network equations are written\n\ | ||
| 16 | in current balance form using rectangular coordinates. It uses the SLEPc\n\ | ||
| 17 | package to calculate the eigenvalues for small signal stability analysis\n\n"; | ||
| 18 | |||
| 19 | /* | ||
| 20 | This example is based on PETSc's ex9bus example (under TS). | ||
| 21 | |||
| 22 | The equations for the stability analysis are described by the DAE | ||
| 23 | |||
| 24 | \dot{x} = f(x,y,t) | ||
| 25 | 0 = g(x,y,t) | ||
| 26 | |||
| 27 | where the generators are described by differential equations, while the algebraic | ||
| 28 | constraints define the network equations. | ||
| 29 | |||
| 30 | The generators are modeled with a 4th order differential equation describing the electrical | ||
| 31 | and mechanical dynamics. Each generator also has an exciter system modeled by 3rd order | ||
| 32 | diff. eqns. describing the exciter, voltage regulator, and the feedback stabilizer | ||
| 33 | mechanism. | ||
| 34 | |||
| 35 | The network equations are described by nodal current balance equations. | ||
| 36 | I(x,y) - Y*V = 0 | ||
| 37 | |||
| 38 | where: | ||
| 39 | I(x,y) is the current injected from generators and loads. | ||
| 40 | Y is the admittance matrix, and | ||
| 41 | V is the voltage vector | ||
| 42 | |||
| 43 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 44 | |||
| 45 | The linearized equations for the eigenvalue analysis are | ||
| 46 | |||
| 47 | \dot{\delta{x}} = f_x*\delta{x} + f_y*\delta{y} | ||
| 48 | 0 = g_x*\delta{x} + g_y*\delta{y} | ||
| 49 | |||
| 50 | This gives the linearized sensitivity matrix | ||
| 51 | A = | f_x f_y | | ||
| 52 | | g_x g_y | | ||
| 53 | |||
| 54 | We are interested in the eigenvalues of the Schur complement of A | ||
| 55 | \hat{A} = f_x - g_x*inv(g_y)*f_y | ||
| 56 | |||
| 57 | Example contributed by: Shrirang Abhyankar | ||
| 58 | */ | ||
| 59 | |||
| 60 | #include <petscdm.h> | ||
| 61 | #include <petscdmda.h> | ||
| 62 | #include <petscdmcomposite.h> | ||
| 63 | #include <slepceps.h> | ||
| 64 | |||
| 65 | #define freq 60 | ||
| 66 | #define w_s (2*PETSC_PI*freq) | ||
| 67 | |||
| 68 | /* Sizes and indices */ | ||
| 69 | const PetscInt nbus = 9; /* Number of network buses */ | ||
| 70 | const PetscInt ngen = 3; /* Number of generators */ | ||
| 71 | const PetscInt nload = 3; /* Number of loads */ | ||
| 72 | const PetscInt gbus[3] = {0,1,2}; /* Buses at which generators are incident */ | ||
| 73 | const PetscInt lbus[3] = {4,5,7}; /* Buses at which loads are incident */ | ||
| 74 | |||
| 75 | /* Generator real and reactive powers (found via loadflow) */ | ||
| 76 | const PetscScalar PG[3] = {0.716786142395021,1.630000000000000,0.850000000000000}; | ||
| 77 | const PetscScalar QG[3] = {0.270702180178785,0.066120127797275,-0.108402221791588}; | ||
| 78 | /* Generator constants */ | ||
| 79 | const PetscScalar H[3] = {23.64,6.4,3.01}; /* Inertia constant */ | ||
| 80 | const PetscScalar Rs[3] = {0.0,0.0,0.0}; /* Stator Resistance */ | ||
| 81 | const PetscScalar Xd[3] = {0.146,0.8958,1.3125}; /* d-axis reactance */ | ||
| 82 | const PetscScalar Xdp[3] = {0.0608,0.1198,0.1813}; /* d-axis transient reactance */ | ||
| 83 | const PetscScalar Xq[3] = {0.0969,0.8645,1.2578}; /* q-axis reactance Xq(1) set to 0.4360, value given in text 0.0969 */ | ||
| 84 | const PetscScalar Xqp[3] = {0.0969,0.1969,0.25}; /* q-axis transient reactance */ | ||
| 85 | const PetscScalar Td0p[3] = {8.96,6.0,5.89}; /* d-axis open circuit time constant */ | ||
| 86 | const PetscScalar Tq0p[3] = {0.31,0.535,0.6}; /* q-axis open circuit time constant */ | ||
| 87 | PetscScalar M[3]; /* M = 2*H/w_s */ | ||
| 88 | PetscScalar D[3]; /* D = 0.1*M */ | ||
| 89 | |||
| 90 | PetscScalar TM[3]; /* Mechanical Torque */ | ||
| 91 | /* Exciter system constants */ | ||
| 92 | const PetscScalar KA[3] = {20.0,20.0,20.0}; /* Voltage regulartor gain constant */ | ||
| 93 | const PetscScalar TA[3] = {0.2,0.2,0.2}; /* Voltage regulator time constant */ | ||
| 94 | const PetscScalar KE[3] = {1.0,1.0,1.0}; /* Exciter gain constant */ | ||
| 95 | const PetscScalar TE[3] = {0.314,0.314,0.314}; /* Exciter time constant */ | ||
| 96 | const PetscScalar KF[3] = {0.063,0.063,0.063}; /* Feedback stabilizer gain constant */ | ||
| 97 | const PetscScalar TF[3] = {0.35,0.35,0.35}; /* Feedback stabilizer time constant */ | ||
| 98 | const PetscScalar k1[3] = {0.0039,0.0039,0.0039}; | ||
| 99 | const PetscScalar k2[3] = {1.555,1.555,1.555}; /* k1 and k2 for calculating the saturation function SE = k1*exp(k2*Efd) */ | ||
| 100 | |||
| 101 | PetscScalar Vref[3]; | ||
| 102 | /* Load constants | ||
| 103 | We use a composite load model that describes the load and reactive powers at each time instant as follows | ||
| 104 | P(t) = \sum\limits_{i=0}^ld_nsegsp \ld_alphap_i*P_D0(\frac{V_m(t)}{V_m0})^\ld_betap_i | ||
| 105 | Q(t) = \sum\limits_{i=0}^ld_nsegsq \ld_alphaq_i*Q_D0(\frac{V_m(t)}{V_m0})^\ld_betaq_i | ||
| 106 | where | ||
| 107 | ld_nsegsp,ld_nsegsq - Number of individual load models for real and reactive power loads | ||
| 108 | ld_alphap,ld_alphap - Percentage contribution (weights) or loads | ||
| 109 | P_D0 - Real power load | ||
| 110 | Q_D0 - Reactive power load | ||
| 111 | V_m(t) - Voltage magnitude at time t | ||
| 112 | V_m0 - Voltage magnitude at t = 0 | ||
| 113 | ld_betap, ld_betaq - exponents describing the load model for real and reactive part | ||
| 114 | |||
| 115 | Note: All loads have the same characteristic currently. | ||
| 116 | */ | ||
| 117 | const PetscScalar PD0[3] = {1.25,0.9,1.0}; | ||
| 118 | const PetscScalar QD0[3] = {0.5,0.3,0.35}; | ||
| 119 | const PetscInt ld_nsegsp[3] = {3,3,3}; | ||
| 120 | const PetscScalar ld_alphap[3] = {0.0,0.0,1.0}; | ||
| 121 | const PetscScalar ld_betap[3] = {2.0,1.0,0.0}; | ||
| 122 | const PetscInt ld_nsegsq[3] = {3,3,3}; | ||
| 123 | const PetscScalar ld_alphaq[3] = {0.0,0.0,1.0}; | ||
| 124 | const PetscScalar ld_betaq[3] = {2.0,1.0,0.0}; | ||
| 125 | |||
| 126 | typedef struct { | ||
| 127 | DM dmgen, dmnet; /* DMs to manage generator and network subsystem */ | ||
| 128 | DM dmpgrid; /* Composite DM to manage the entire power grid */ | ||
| 129 | Mat Ybus; /* Network admittance matrix */ | ||
| 130 | Vec V0; /* Initial voltage vector (Power flow solution) */ | ||
| 131 | PetscInt neqs_gen,neqs_net,neqs_pgrid; | ||
| 132 | IS is_diff; /* indices for differential equations */ | ||
| 133 | IS is_alg; /* indices for algebraic equations */ | ||
| 134 | } Userctx; | ||
| 135 | |||
| 136 | /* Converts from machine frame (dq) to network (phase a real,imag) reference frame */ | ||
| 137 | ✗ | PetscErrorCode dq2ri(PetscScalar Fd,PetscScalar Fq,PetscScalar delta,PetscScalar *Fr,PetscScalar *Fi) | |
| 138 | { | ||
| 139 | ✗ | PetscFunctionBegin; | |
| 140 | ✗ | *Fr = Fd*PetscSinScalar(delta) + Fq*PetscCosScalar(delta); | |
| 141 | ✗ | *Fi = -Fd*PetscCosScalar(delta) + Fq*PetscSinScalar(delta); | |
| 142 | ✗ | PetscFunctionReturn(PETSC_SUCCESS); | |
| 143 | } | ||
| 144 | |||
| 145 | /* Converts from network frame ([phase a real,imag) to machine (dq) reference frame */ | ||
| 146 | 15 | PetscErrorCode ri2dq(PetscScalar Fr,PetscScalar Fi,PetscScalar delta,PetscScalar *Fd,PetscScalar *Fq) | |
| 147 | { | ||
| 148 |
1/2✓ Branch 0 taken 1 times.
✗ Branch 1 not taken.
|
15 | PetscFunctionBegin; |
| 149 | 15 | *Fd = Fr*PetscSinScalar(delta) - Fi*PetscCosScalar(delta); | |
| 150 | 15 | *Fq = Fr*PetscCosScalar(delta) + Fi*PetscSinScalar(delta); | |
| 151 |
6/12✓ Branch 0 taken 1 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 1 times.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 1 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 1 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 1 times.
|
15 | PetscFunctionReturn(PETSC_SUCCESS); |
| 152 | } | ||
| 153 | |||
| 154 | 5 | PetscErrorCode SetInitialGuess(Vec X,Userctx *user) | |
| 155 | { | ||
| 156 | 5 | Vec Xgen,Xnet; | |
| 157 | 5 | PetscScalar *xgen,*xnet; | |
| 158 | 5 | PetscInt i,idx=0; | |
| 159 | 5 | PetscScalar Vr,Vi,IGr,IGi,Vm,Vm2; | |
| 160 | 5 | PetscScalar Eqp,Edp,delta; | |
| 161 | 5 | PetscScalar Efd,RF,VR; /* Exciter variables */ | |
| 162 | 5 | PetscScalar Id,Iq; /* Generator dq axis currents */ | |
| 163 | 5 | PetscScalar theta,Vd,Vq,SE; | |
| 164 | |||
| 165 |
1/2✓ Branch 0 taken 1 times.
✗ Branch 1 not taken.
|
5 | PetscFunctionBegin; |
| 166 | 5 | M[0] = 2*H[0]/w_s; M[1] = 2*H[1]/w_s; M[2] = 2*H[2]/w_s; | |
| 167 | /* D[0] = 0.1*M[0]; D[1] = 0.1*M[1]; D[2] = 0.1*M[2]; | ||
| 168 | */ | ||
| 169 | 5 | D[0] = D[1] = D[2] = 0.0; | |
| 170 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet)); |
| 171 | |||
| 172 | /* Network subsystem initialization */ | ||
| 173 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(VecCopy(user->V0,Xnet)); |
| 174 | |||
| 175 | /* Generator subsystem initialization */ | ||
| 176 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(VecGetArray(Xgen,&xgen)); |
| 177 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(VecGetArray(Xnet,&xnet)); |
| 178 | |||
| 179 |
2/2✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
|
20 | for (i=0; i < ngen; i++) { |
| 180 | 15 | Vr = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */ | |
| 181 | 15 | Vi = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */ | |
| 182 | 15 | Vm = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm; | |
| 183 | 15 | IGr = (Vr*PG[i] + Vi*QG[i])/Vm2; | |
| 184 | 15 | IGi = (Vi*PG[i] - Vr*QG[i])/Vm2; | |
| 185 | |||
| 186 | 15 | delta = atan2(Vi+Xq[i]*IGr,Vr-Xq[i]*IGi); /* Machine angle */ | |
| 187 | |||
| 188 | 15 | theta = PETSC_PI/2.0 - delta; | |
| 189 | |||
| 190 | 15 | Id = IGr*PetscCosScalar(theta) - IGi*PetscSinScalar(theta); /* d-axis stator current */ | |
| 191 | 15 | Iq = IGr*PetscSinScalar(theta) + IGi*PetscCosScalar(theta); /* q-axis stator current */ | |
| 192 | |||
| 193 | 15 | Vd = Vr*PetscCosScalar(theta) - Vi*PetscSinScalar(theta); | |
| 194 | 15 | Vq = Vr*PetscSinScalar(theta) + Vi*PetscCosScalar(theta); | |
| 195 | |||
| 196 | 15 | Edp = Vd + Rs[i]*Id - Xqp[i]*Iq; /* d-axis transient EMF */ | |
| 197 | 15 | Eqp = Vq + Rs[i]*Iq + Xdp[i]*Id; /* q-axis transient EMF */ | |
| 198 | |||
| 199 | 15 | TM[i] = PG[i]; | |
| 200 | |||
| 201 | /* The generator variables are ordered as [Eqp,Edp,delta,w,Id,Iq] */ | ||
| 202 | 15 | xgen[idx] = Eqp; | |
| 203 | 15 | xgen[idx+1] = Edp; | |
| 204 | 15 | xgen[idx+2] = delta; | |
| 205 | 15 | xgen[idx+3] = w_s; | |
| 206 | |||
| 207 | 15 | idx = idx + 4; | |
| 208 | |||
| 209 | 15 | xgen[idx] = Id; | |
| 210 | 15 | xgen[idx+1] = Iq; | |
| 211 | |||
| 212 | 15 | idx = idx + 2; | |
| 213 | |||
| 214 | /* Exciter */ | ||
| 215 | 15 | Efd = Eqp + (Xd[i] - Xdp[i])*Id; | |
| 216 | 15 | SE = k1[i]*PetscExpScalar(k2[i]*Efd); | |
| 217 | 15 | VR = KE[i]*Efd + SE; | |
| 218 | 15 | RF = KF[i]*Efd/TF[i]; | |
| 219 | |||
| 220 | 15 | xgen[idx] = Efd; | |
| 221 | 15 | xgen[idx+1] = RF; | |
| 222 | 15 | xgen[idx+2] = VR; | |
| 223 | |||
| 224 | 15 | Vref[i] = Vm + (VR/KA[i]); | |
| 225 | |||
| 226 | 15 | idx = idx + 3; | |
| 227 | } | ||
| 228 | |||
| 229 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(VecRestoreArray(Xgen,&xgen)); |
| 230 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(VecRestoreArray(Xnet,&xnet)); |
| 231 | |||
| 232 | /* PetscCall(VecView(Xgen,0)); */ | ||
| 233 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(DMCompositeGather(user->dmpgrid,INSERT_VALUES,X,Xgen,Xnet)); |
| 234 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet)); |
| 235 |
6/12✓ Branch 0 taken 1 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 1 times.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 1 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 1 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 1 times.
|
1 | PetscFunctionReturn(PETSC_SUCCESS); |
| 236 | } | ||
| 237 | |||
| 238 | 5 | PetscErrorCode PreallocateJacobian(Mat J,Userctx *user) | |
| 239 | { | ||
| 240 | 5 | PetscInt *d_nnz; | |
| 241 | 5 | PetscInt i,idx=0,start=0; | |
| 242 | 5 | PetscInt ncols; | |
| 243 | |||
| 244 |
1/2✓ Branch 0 taken 1 times.
✗ Branch 1 not taken.
|
5 | PetscFunctionBegin; |
| 245 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(PetscMalloc1(user->neqs_pgrid,&d_nnz)); |
| 246 |
2/2✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
|
230 | for (i=0; i<user->neqs_pgrid; i++) d_nnz[i] = 0; |
| 247 | /* Generator subsystem */ | ||
| 248 |
2/2✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
|
20 | for (i=0; i < ngen; i++) { |
| 249 | |||
| 250 | 15 | d_nnz[idx] += 3; | |
| 251 | 15 | d_nnz[idx+1] += 2; | |
| 252 | 15 | d_nnz[idx+2] += 2; | |
| 253 | 15 | d_nnz[idx+3] += 5; | |
| 254 | 15 | d_nnz[idx+4] += 6; | |
| 255 | 15 | d_nnz[idx+5] += 6; | |
| 256 | |||
| 257 | 15 | d_nnz[user->neqs_gen+2*gbus[i]] += 3; | |
| 258 | 15 | d_nnz[user->neqs_gen+2*gbus[i]+1] += 3; | |
| 259 | |||
| 260 | 15 | d_nnz[idx+6] += 2; | |
| 261 | 15 | d_nnz[idx+7] += 2; | |
| 262 | 15 | d_nnz[idx+8] += 5; | |
| 263 | |||
| 264 | 15 | idx = idx + 9; | |
| 265 | } | ||
| 266 | |||
| 267 | 5 | start = user->neqs_gen; | |
| 268 | |||
| 269 |
2/2✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
|
50 | for (i=0; i < nbus; i++) { |
| 270 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
45 | PetscCall(MatGetRow(user->Ybus,2*i,&ncols,NULL,NULL)); |
| 271 | 45 | d_nnz[start+2*i] += ncols; | |
| 272 | 45 | d_nnz[start+2*i+1] += ncols; | |
| 273 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
45 | PetscCall(MatRestoreRow(user->Ybus,2*i,&ncols,NULL,NULL)); |
| 274 | } | ||
| 275 | |||
| 276 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(MatSeqAIJSetPreallocation(J,0,d_nnz)); |
| 277 | |||
| 278 |
5/8✓ Branch 0 taken 5 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 1 times.
✓ Branch 3 taken 4 times.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 1 times.
|
5 | PetscCall(PetscFree(d_nnz)); |
| 279 |
6/12✓ Branch 0 taken 1 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 1 times.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 1 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 1 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 1 times.
|
1 | PetscFunctionReturn(PETSC_SUCCESS); |
| 280 | } | ||
| 281 | |||
| 282 | /* | ||
| 283 | J = [-df_dx, -df_dy | ||
| 284 | dg_dx, dg_dy] | ||
| 285 | */ | ||
| 286 | 5 | PetscErrorCode ResidualJacobian(Vec X,Mat J,void *ctx) | |
| 287 | { | ||
| 288 | 5 | Userctx *user=(Userctx*)ctx; | |
| 289 | 5 | Vec Xgen,Xnet; | |
| 290 | 5 | PetscScalar *xgen,*xnet; | |
| 291 | 5 | PetscInt i,idx=0; | |
| 292 | 5 | PetscScalar Vr,Vi,Vm,Vm2; | |
| 293 | 5 | PetscScalar Eqp,Edp,delta; /* Generator variables */ | |
| 294 | 5 | PetscScalar Efd; | |
| 295 | 5 | PetscScalar Id,Iq; /* Generator dq axis currents */ | |
| 296 | 5 | PetscScalar Vd,Vq; | |
| 297 | 5 | PetscScalar val[10]; | |
| 298 | 5 | PetscInt row[2],col[10]; | |
| 299 | 5 | PetscInt net_start=user->neqs_gen; | |
| 300 | 5 | PetscScalar Zdq_inv[4],det; | |
| 301 | 5 | PetscScalar dVd_dVr,dVd_dVi,dVq_dVr,dVq_dVi,dVd_ddelta,dVq_ddelta; | |
| 302 | 5 | PetscScalar dIGr_ddelta,dIGi_ddelta,dIGr_dId,dIGr_dIq,dIGi_dId,dIGi_dIq; | |
| 303 | 5 | PetscScalar dSE_dEfd; | |
| 304 | 5 | PetscScalar dVm_dVd,dVm_dVq,dVm_dVr,dVm_dVi; | |
| 305 | 5 | PetscInt ncols; | |
| 306 | 5 | const PetscInt *cols; | |
| 307 | 5 | const PetscScalar *yvals; | |
| 308 | 5 | PetscInt k; | |
| 309 | 5 | PetscScalar PD,QD,Vm0,*v0,Vm4; | |
| 310 | 5 | PetscScalar dPD_dVr,dPD_dVi,dQD_dVr,dQD_dVi; | |
| 311 | 5 | PetscScalar dIDr_dVr,dIDr_dVi,dIDi_dVr,dIDi_dVi; | |
| 312 | |||
| 313 |
1/2✓ Branch 0 taken 1 times.
✗ Branch 1 not taken.
|
5 | PetscFunctionBegin; |
| 314 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(MatZeroEntries(J)); |
| 315 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet)); |
| 316 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(DMCompositeScatter(user->dmpgrid,X,Xgen,Xnet)); |
| 317 | |||
| 318 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(VecGetArray(Xgen,&xgen)); |
| 319 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(VecGetArray(Xnet,&xnet)); |
| 320 | |||
| 321 | /* Generator subsystem */ | ||
| 322 |
2/2✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
|
20 | for (i=0; i < ngen; i++) { |
| 323 | 15 | Eqp = xgen[idx]; | |
| 324 | 15 | Edp = xgen[idx+1]; | |
| 325 | 15 | delta = xgen[idx+2]; | |
| 326 | 15 | Id = xgen[idx+4]; | |
| 327 | 15 | Iq = xgen[idx+5]; | |
| 328 | 15 | Efd = xgen[idx+6]; | |
| 329 | |||
| 330 | /* fgen[idx] = (Eqp + (Xd[i] - Xdp[i])*Id - Efd)/Td0p[i]; */ | ||
| 331 | 15 | row[0] = idx; | |
| 332 | 15 | col[0] = idx; col[1] = idx+4; col[2] = idx+6; | |
| 333 | 15 | val[0] = 1/ Td0p[i]; val[1] = (Xd[i] - Xdp[i])/ Td0p[i]; val[2] = -1/Td0p[i]; | |
| 334 | |||
| 335 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
15 | PetscCall(MatSetValues(J,1,row,3,col,val,INSERT_VALUES)); |
| 336 | |||
| 337 | /* fgen[idx+1] = (Edp - (Xq[i] - Xqp[i])*Iq)/Tq0p[i]; */ | ||
| 338 | 15 | row[0] = idx + 1; | |
| 339 | 15 | col[0] = idx + 1; col[1] = idx+5; | |
| 340 | 15 | val[0] = 1/Tq0p[i]; val[1] = -(Xq[i] - Xqp[i])/Tq0p[i]; | |
| 341 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
15 | PetscCall(MatSetValues(J,1,row,2,col,val,INSERT_VALUES)); |
| 342 | |||
| 343 | /* fgen[idx+2] = - w + w_s; */ | ||
| 344 | 15 | row[0] = idx + 2; | |
| 345 | 15 | col[0] = idx + 2; col[1] = idx + 3; | |
| 346 | 15 | val[0] = 0; val[1] = -1; | |
| 347 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
15 | PetscCall(MatSetValues(J,1,row,2,col,val,INSERT_VALUES)); |
| 348 | |||
| 349 | /* fgen[idx+3] = (-TM[i] + Edp*Id + Eqp*Iq + (Xqp[i] - Xdp[i])*Id*Iq + D[i]*(w - w_s))/M[i]; */ | ||
| 350 | 15 | row[0] = idx + 3; | |
| 351 | 15 | col[0] = idx; col[1] = idx + 1; col[2] = idx + 3; col[3] = idx + 4; col[4] = idx + 5; | |
| 352 | 15 | 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]; | |
| 353 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
15 | PetscCall(MatSetValues(J,1,row,5,col,val,INSERT_VALUES)); |
| 354 | |||
| 355 | 15 | Vr = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */ | |
| 356 | 15 | Vi = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */ | |
| 357 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
15 | PetscCall(ri2dq(Vr,Vi,delta,&Vd,&Vq)); |
| 358 | |||
| 359 | 15 | det = Rs[i]*Rs[i] + Xdp[i]*Xqp[i]; | |
| 360 | |||
| 361 | 15 | Zdq_inv[0] = Rs[i]/det; | |
| 362 | 15 | Zdq_inv[1] = Xqp[i]/det; | |
| 363 | 15 | Zdq_inv[2] = -Xdp[i]/det; | |
| 364 | 15 | Zdq_inv[3] = Rs[i]/det; | |
| 365 | |||
| 366 | 15 | dVd_dVr = PetscSinScalar(delta); dVd_dVi = -PetscCosScalar(delta); | |
| 367 | 15 | dVq_dVr = PetscCosScalar(delta); dVq_dVi = PetscSinScalar(delta); | |
| 368 | 15 | dVd_ddelta = Vr*PetscCosScalar(delta) + Vi*PetscSinScalar(delta); | |
| 369 | 15 | dVq_ddelta = -Vr*PetscSinScalar(delta) + Vi*PetscCosScalar(delta); | |
| 370 | |||
| 371 | /* fgen[idx+4] = Zdq_inv[0]*(-Edp + Vd) + Zdq_inv[1]*(-Eqp + Vq) + Id; */ | ||
| 372 | 15 | row[0] = idx+4; | |
| 373 | 15 | col[0] = idx; col[1] = idx+1; col[2] = idx + 2; | |
| 374 | 15 | val[0] = -Zdq_inv[1]; val[1] = -Zdq_inv[0]; val[2] = Zdq_inv[0]*dVd_ddelta + Zdq_inv[1]*dVq_ddelta; | |
| 375 | 15 | col[3] = idx + 4; col[4] = net_start+2*gbus[i]; col[5] = net_start + 2*gbus[i]+1; | |
| 376 | 15 | 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; | |
| 377 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
15 | PetscCall(MatSetValues(J,1,row,6,col,val,INSERT_VALUES)); |
| 378 | |||
| 379 | /* fgen[idx+5] = Zdq_inv[2]*(-Edp + Vd) + Zdq_inv[3]*(-Eqp + Vq) + Iq; */ | ||
| 380 | 15 | row[0] = idx+5; | |
| 381 | 15 | col[0] = idx; col[1] = idx+1; col[2] = idx + 2; | |
| 382 | 15 | val[0] = -Zdq_inv[3]; val[1] = -Zdq_inv[2]; val[2] = Zdq_inv[2]*dVd_ddelta + Zdq_inv[3]*dVq_ddelta; | |
| 383 | 15 | col[3] = idx + 5; col[4] = net_start+2*gbus[i]; col[5] = net_start + 2*gbus[i]+1; | |
| 384 | 15 | 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; | |
| 385 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
15 | PetscCall(MatSetValues(J,1,row,6,col,val,INSERT_VALUES)); |
| 386 | |||
| 387 | 15 | dIGr_ddelta = Id*PetscCosScalar(delta) - Iq*PetscSinScalar(delta); | |
| 388 | 15 | dIGi_ddelta = Id*PetscSinScalar(delta) + Iq*PetscCosScalar(delta); | |
| 389 | 15 | dIGr_dId = PetscSinScalar(delta); dIGr_dIq = PetscCosScalar(delta); | |
| 390 | 15 | dIGi_dId = -PetscCosScalar(delta); dIGi_dIq = PetscSinScalar(delta); | |
| 391 | |||
| 392 | /* fnet[2*gbus[i]] -= IGi; */ | ||
| 393 | 15 | row[0] = net_start + 2*gbus[i]; | |
| 394 | 15 | col[0] = idx+2; col[1] = idx + 4; col[2] = idx + 5; | |
| 395 | 15 | val[0] = -dIGi_ddelta; val[1] = -dIGi_dId; val[2] = -dIGi_dIq; | |
| 396 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
15 | PetscCall(MatSetValues(J,1,row,3,col,val,INSERT_VALUES)); |
| 397 | |||
| 398 | /* fnet[2*gbus[i]+1] -= IGr; */ | ||
| 399 | 15 | row[0] = net_start + 2*gbus[i]+1; | |
| 400 | 15 | col[0] = idx+2; col[1] = idx + 4; col[2] = idx + 5; | |
| 401 | 15 | val[0] = -dIGr_ddelta; val[1] = -dIGr_dId; val[2] = -dIGr_dIq; | |
| 402 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
15 | PetscCall(MatSetValues(J,1,row,3,col,val,INSERT_VALUES)); |
| 403 | |||
| 404 | 15 | Vm = PetscSqrtScalar(Vd*Vd + Vq*Vq); Vm2 = Vm*Vm; | |
| 405 | |||
| 406 | /* fgen[idx+6] = (KE[i]*Efd + SE - VR)/TE[i]; */ | ||
| 407 | /* SE = k1[i]*PetscExpScalar(k2[i]*Efd); */ | ||
| 408 | |||
| 409 | 15 | dSE_dEfd = k1[i]*k2[i]*PetscExpScalar(k2[i]*Efd); | |
| 410 | |||
| 411 | 15 | row[0] = idx + 6; | |
| 412 | 15 | col[0] = idx + 6; col[1] = idx + 8; | |
| 413 | 15 | val[0] = (KE[i] + dSE_dEfd)/TE[i]; val[1] = -1/TE[i]; | |
| 414 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
15 | PetscCall(MatSetValues(J,1,row,2,col,val,INSERT_VALUES)); |
| 415 | |||
| 416 | /* Exciter differential equations */ | ||
| 417 | |||
| 418 | /* fgen[idx+7] = (RF - KF[i]*Efd/TF[i])/TF[i]; */ | ||
| 419 | 15 | row[0] = idx + 7; | |
| 420 | 15 | col[0] = idx + 6; col[1] = idx + 7; | |
| 421 | 15 | val[0] = (-KF[i]/TF[i])/TF[i]; val[1] = 1/TF[i]; | |
| 422 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
15 | PetscCall(MatSetValues(J,1,row,2,col,val,INSERT_VALUES)); |
| 423 | |||
| 424 | /* fgen[idx+8] = (VR - KA[i]*RF + KA[i]*KF[i]*Efd/TF[i] - KA[i]*(Vref[i] - Vm))/TA[i]; */ | ||
| 425 | /* Vm = (Vd^2 + Vq^2)^0.5; */ | ||
| 426 | |||
| 427 | 15 | dVm_dVd = Vd/Vm; dVm_dVq = Vq/Vm; | |
| 428 | 15 | dVm_dVr = dVm_dVd*dVd_dVr + dVm_dVq*dVq_dVr; | |
| 429 | 15 | dVm_dVi = dVm_dVd*dVd_dVi + dVm_dVq*dVq_dVi; | |
| 430 | 15 | row[0] = idx + 8; | |
| 431 | 15 | col[0] = idx + 6; col[1] = idx + 7; col[2] = idx + 8; | |
| 432 | 15 | val[0] = (KA[i]*KF[i]/TF[i])/TA[i]; val[1] = -KA[i]/TA[i]; val[2] = 1/TA[i]; | |
| 433 | 15 | col[3] = net_start + 2*gbus[i]; col[4] = net_start + 2*gbus[i]+1; | |
| 434 | 15 | val[3] = KA[i]*dVm_dVr/TA[i]; val[4] = KA[i]*dVm_dVi/TA[i]; | |
| 435 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
15 | PetscCall(MatSetValues(J,1,row,5,col,val,INSERT_VALUES)); |
| 436 | 15 | idx = idx + 9; | |
| 437 | } | ||
| 438 | |||
| 439 |
2/2✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
|
50 | for (i=0; i<nbus; i++) { |
| 440 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
45 | PetscCall(MatGetRow(user->Ybus,2*i,&ncols,&cols,&yvals)); |
| 441 | 45 | row[0] = net_start + 2*i; | |
| 442 |
2/2✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
|
315 | for (k=0; k<ncols; k++) { |
| 443 | 270 | col[k] = net_start + cols[k]; | |
| 444 | 270 | val[k] = yvals[k]; | |
| 445 | } | ||
| 446 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
45 | PetscCall(MatSetValues(J,1,row,ncols,col,val,INSERT_VALUES)); |
| 447 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
45 | PetscCall(MatRestoreRow(user->Ybus,2*i,&ncols,&cols,&yvals)); |
| 448 | |||
| 449 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
45 | PetscCall(MatGetRow(user->Ybus,2*i+1,&ncols,&cols,&yvals)); |
| 450 | 45 | row[0] = net_start + 2*i+1; | |
| 451 |
2/2✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
|
315 | for (k=0; k<ncols; k++) { |
| 452 | 270 | col[k] = net_start + cols[k]; | |
| 453 | 270 | val[k] = yvals[k]; | |
| 454 | } | ||
| 455 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
45 | PetscCall(MatSetValues(J,1,row,ncols,col,val,INSERT_VALUES)); |
| 456 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
45 | PetscCall(MatRestoreRow(user->Ybus,2*i+1,&ncols,&cols,&yvals)); |
| 457 | } | ||
| 458 | |||
| 459 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(MatAssemblyBegin(J,MAT_FLUSH_ASSEMBLY)); |
| 460 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(MatAssemblyEnd(J,MAT_FLUSH_ASSEMBLY)); |
| 461 | |||
| 462 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(VecGetArray(user->V0,&v0)); |
| 463 |
2/2✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
|
20 | for (i=0; i < nload; i++) { |
| 464 | 15 | Vr = xnet[2*lbus[i]]; /* Real part of load bus voltage */ | |
| 465 | 15 | Vi = xnet[2*lbus[i]+1]; /* Imaginary part of the load bus voltage */ | |
| 466 | 15 | Vm = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm; Vm4 = Vm2*Vm2; | |
| 467 | 15 | Vm0 = PetscSqrtScalar(v0[2*lbus[i]]*v0[2*lbus[i]] + v0[2*lbus[i]+1]*v0[2*lbus[i]+1]); | |
| 468 | 15 | PD = QD = 0.0; | |
| 469 | 15 | dPD_dVr = dPD_dVi = dQD_dVr = dQD_dVi = 0.0; | |
| 470 |
2/2✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
|
60 | for (k=0; k < ld_nsegsp[i]; k++) { |
| 471 | 45 | PD += ld_alphap[k]*PD0[i]*PetscPowScalar((Vm/Vm0),ld_betap[k]); | |
| 472 | 45 | dPD_dVr += ld_alphap[k]*ld_betap[k]*PD0[i]*PetscPowScalar((1/Vm0),ld_betap[k])*Vr*PetscPowScalar(Vm,(ld_betap[k]-2)); | |
| 473 | 45 | dPD_dVi += ld_alphap[k]*ld_betap[k]*PD0[i]*PetscPowScalar((1/Vm0),ld_betap[k])*Vi*PetscPowScalar(Vm,(ld_betap[k]-2)); | |
| 474 | } | ||
| 475 |
2/2✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
|
60 | for (k=0; k < ld_nsegsq[i]; k++) { |
| 476 | 45 | QD += ld_alphaq[k]*QD0[i]*PetscPowScalar((Vm/Vm0),ld_betaq[k]); | |
| 477 | 45 | dQD_dVr += ld_alphaq[k]*ld_betaq[k]*QD0[i]*PetscPowScalar((1/Vm0),ld_betaq[k])*Vr*PetscPowScalar(Vm,(ld_betaq[k]-2)); | |
| 478 | 45 | dQD_dVi += ld_alphaq[k]*ld_betaq[k]*QD0[i]*PetscPowScalar((1/Vm0),ld_betaq[k])*Vi*PetscPowScalar(Vm,(ld_betaq[k]-2)); | |
| 479 | } | ||
| 480 | |||
| 481 | /* IDr = (PD*Vr + QD*Vi)/Vm2; */ | ||
| 482 | /* IDi = (-QD*Vr + PD*Vi)/Vm2; */ | ||
| 483 | |||
| 484 | 15 | dIDr_dVr = (dPD_dVr*Vr + dQD_dVr*Vi + PD)/Vm2 - ((PD*Vr + QD*Vi)*2*Vr)/Vm4; | |
| 485 | 15 | dIDr_dVi = (dPD_dVi*Vr + dQD_dVi*Vi + QD)/Vm2 - ((PD*Vr + QD*Vi)*2*Vi)/Vm4; | |
| 486 | |||
| 487 | 15 | dIDi_dVr = (-dQD_dVr*Vr + dPD_dVr*Vi - QD)/Vm2 - ((-QD*Vr + PD*Vi)*2*Vr)/Vm4; | |
| 488 | 15 | dIDi_dVi = (-dQD_dVi*Vr + dPD_dVi*Vi + PD)/Vm2 - ((-QD*Vr + PD*Vi)*2*Vi)/Vm4; | |
| 489 | |||
| 490 | /* fnet[2*lbus[i]] += IDi; */ | ||
| 491 | 15 | row[0] = net_start + 2*lbus[i]; | |
| 492 | 15 | col[0] = net_start + 2*lbus[i]; col[1] = net_start + 2*lbus[i]+1; | |
| 493 | 15 | val[0] = dIDi_dVr; val[1] = dIDi_dVi; | |
| 494 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
15 | PetscCall(MatSetValues(J,1,row,2,col,val,ADD_VALUES)); |
| 495 | /* fnet[2*lbus[i]+1] += IDr; */ | ||
| 496 | 15 | row[0] = net_start + 2*lbus[i]+1; | |
| 497 | 15 | col[0] = net_start + 2*lbus[i]; col[1] = net_start + 2*lbus[i]+1; | |
| 498 | 15 | val[0] = dIDr_dVr; val[1] = dIDr_dVi; | |
| 499 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
15 | PetscCall(MatSetValues(J,1,row,2,col,val,ADD_VALUES)); |
| 500 | } | ||
| 501 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(VecRestoreArray(user->V0,&v0)); |
| 502 | |||
| 503 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(VecRestoreArray(Xgen,&xgen)); |
| 504 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(VecRestoreArray(Xnet,&xnet)); |
| 505 | |||
| 506 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet)); |
| 507 | |||
| 508 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY)); |
| 509 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY)); |
| 510 |
6/12✓ Branch 0 taken 1 times.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✓ Branch 3 taken 1 times.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 1 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 1 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 1 times.
|
1 | PetscFunctionReturn(PETSC_SUCCESS); |
| 511 | } | ||
| 512 | |||
| 513 | 5 | int main(int argc,char **argv) | |
| 514 | { | ||
| 515 | 5 | EPS eps; | |
| 516 | 5 | EPSType type; | |
| 517 | 5 | PetscMPIInt size; | |
| 518 | 5 | Userctx user; | |
| 519 | 5 | PetscViewer Xview,Ybusview; | |
| 520 | 5 | Vec X,Xr,Xi; | |
| 521 | 5 | Mat J,Jred=NULL; | |
| 522 | 5 | IS is0,is1; | |
| 523 | 5 | PetscInt i,*idx2,its,nev,nconv; | |
| 524 | 5 | PetscReal error,re,im; | |
| 525 | 5 | PetscScalar kr,ki; | |
| 526 | 5 | PetscBool terse; | |
| 527 | |||
| 528 |
1/2✓ Branch 0 taken 1 times.
✗ Branch 1 not taken.
|
5 | PetscFunctionBeginUser; |
| 529 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(SlepcInitialize(&argc,&argv,NULL,help)); |
| 530 |
14/28✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✓ Branch 3 taken 4 times.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 1 times.
✗ Branch 7 not taken.
✓ Branch 8 taken 1 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 1 times.
✓ Branch 12 taken 1 times.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 1 times.
✓ Branch 16 taken 1 times.
✗ Branch 17 not taken.
✓ Branch 18 taken 1 times.
✗ Branch 19 not taken.
✗ Branch 20 not taken.
✓ Branch 21 taken 1 times.
✗ Branch 22 not taken.
✗ Branch 23 not taken.
✗ Branch 24 not taken.
✓ Branch 25 taken 1 times.
✗ Branch 26 not taken.
✗ Branch 27 not taken.
|
5 | PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); |
| 531 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 5 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
5 | PetscCheck(size==1,PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs"); |
| 532 | /* show detailed info unless -terse option is given by user */ | ||
| 533 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse)); |
| 534 | |||
| 535 | 5 | user.neqs_gen = 9*ngen; /* # eqs. for generator subsystem */ | |
| 536 | 5 | user.neqs_net = 2*nbus; /* # eqs. for network subsystem */ | |
| 537 | 5 | user.neqs_pgrid = user.neqs_gen + user.neqs_net; | |
| 538 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nStability analysis in a network with %" PetscInt_FMT " buses and %" PetscInt_FMT " generators\n\n",nbus,ngen)); |
| 539 | |||
| 540 | /* Create indices for differential and algebraic equations */ | ||
| 541 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(PetscMalloc1(7*ngen,&idx2)); |
| 542 |
2/2✓ Branch 0 taken 5 times.
✓ Branch 1 taken 5 times.
|
20 | for (i=0; i<ngen; i++) { |
| 543 | 15 | 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; | |
| 544 | 15 | idx2[7*i+4] = 9*i+6; idx2[7*i+5] = 9*i+7; idx2[7*i+6] = 9*i+8; | |
| 545 | } | ||
| 546 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,7*ngen,idx2,PETSC_COPY_VALUES,&user.is_diff)); |
| 547 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(ISComplement(user.is_diff,0,user.neqs_pgrid,&user.is_alg)); |
| 548 |
5/8✓ Branch 0 taken 5 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 1 times.
✓ Branch 3 taken 4 times.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 1 times.
|
5 | PetscCall(PetscFree(idx2)); |
| 549 | |||
| 550 | /* Read initial voltage vector and Ybus */ | ||
| 551 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"X.bin",FILE_MODE_READ,&Xview)); |
| 552 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"Ybus.bin",FILE_MODE_READ,&Ybusview)); |
| 553 | |||
| 554 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(VecCreate(PETSC_COMM_WORLD,&user.V0)); |
| 555 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(VecSetSizes(user.V0,PETSC_DECIDE,user.neqs_net)); |
| 556 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(VecLoad(user.V0,Xview)); |
| 557 | |||
| 558 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(MatCreate(PETSC_COMM_WORLD,&user.Ybus)); |
| 559 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(MatSetSizes(user.Ybus,PETSC_DECIDE,PETSC_DECIDE,user.neqs_net,user.neqs_net)); |
| 560 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(MatSetType(user.Ybus,MATBAIJ)); |
| 561 | /* PetscCall(MatSetBlockSize(user.Ybus,2)); */ | ||
| 562 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(MatLoad(user.Ybus,Ybusview)); |
| 563 | |||
| 564 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(PetscViewerDestroy(&Xview)); |
| 565 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(PetscViewerDestroy(&Ybusview)); |
| 566 | |||
| 567 | /* Create DMs for generator and network subsystems */ | ||
| 568 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,user.neqs_gen,1,1,NULL,&user.dmgen)); |
| 569 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(DMSetOptionsPrefix(user.dmgen,"dmgen_")); |
| 570 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(DMSetFromOptions(user.dmgen)); |
| 571 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(DMSetUp(user.dmgen)); |
| 572 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,user.neqs_net,1,1,NULL,&user.dmnet)); |
| 573 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(DMSetOptionsPrefix(user.dmnet,"dmnet_")); |
| 574 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(DMSetFromOptions(user.dmnet)); |
| 575 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(DMSetUp(user.dmnet)); |
| 576 | |||
| 577 | /* Create a composite DM packer and add the two DMs */ | ||
| 578 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(DMCompositeCreate(PETSC_COMM_WORLD,&user.dmpgrid)); |
| 579 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(DMSetOptionsPrefix(user.dmpgrid,"pgrid_")); |
| 580 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(DMCompositeAddDM(user.dmpgrid,user.dmgen)); |
| 581 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(DMCompositeAddDM(user.dmpgrid,user.dmnet)); |
| 582 | |||
| 583 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(DMCreateGlobalVector(user.dmpgrid,&X)); |
| 584 | |||
| 585 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(MatCreate(PETSC_COMM_WORLD,&J)); |
| 586 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,user.neqs_pgrid,user.neqs_pgrid)); |
| 587 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(MatSetFromOptions(J)); |
| 588 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(PreallocateJacobian(J,&user)); |
| 589 | |||
| 590 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 591 | Set initial conditions | ||
| 592 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
| 593 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(SetInitialGuess(X,&user)); |
| 594 | |||
| 595 | /* Form Jacobian */ | ||
| 596 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(ResidualJacobian(X,J,(void*)&user)); |
| 597 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(MatScale(J,-1)); |
| 598 | 5 | is0 = user.is_diff; | |
| 599 | 5 | is1 = user.is_alg; | |
| 600 | |||
| 601 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(MatGetSchurComplement(J,is1,is1,is0,is0,MAT_IGNORE_MATRIX,NULL,MAT_SCHUR_COMPLEMENT_AINV_DIAG,MAT_INITIAL_MATRIX,&Jred)); |
| 602 | |||
| 603 |
1/8✗ Branch 0 not taken.
✓ Branch 1 taken 5 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
5 | if (!terse) PetscCall(MatView(Jred,NULL)); |
| 604 | |||
| 605 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(MatCreateVecs(Jred,NULL,&Xr)); |
| 606 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(MatCreateVecs(Jred,NULL,&Xi)); |
| 607 | |||
| 608 | /* Create the eigensolver and set the various options */ | ||
| 609 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps)); |
| 610 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(EPSSetOperators(eps,Jred,NULL)); |
| 611 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(EPSSetProblemType(eps,EPS_NHEP)); |
| 612 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(EPSSetFromOptions(eps)); |
| 613 | |||
| 614 | /* Solve the eigenvalue problem */ | ||
| 615 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(EPSSolve(eps)); |
| 616 | |||
| 617 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(EPSGetIterationNumber(eps,&its)); |
| 618 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the eigensolver: %" PetscInt_FMT "\n",its)); |
| 619 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(EPSGetType(eps,&type)); |
| 620 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n", type)); |
| 621 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(EPSGetDimensions(eps,&nev,NULL,NULL)); |
| 622 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev)); |
| 623 | |||
| 624 | /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - | ||
| 625 | Display solution and clean up | ||
| 626 | - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ | ||
| 627 |
5/8✓ Branch 0 taken 5 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 1 times.
✓ Branch 3 taken 4 times.
✓ Branch 4 taken 1 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 7 taken 1 times.
|
5 | if (terse) PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL)); |
| 628 | else { | ||
| 629 | /* Get number of converged approximate eigenpairs */ | ||
| 630 | ✗ | PetscCall(EPSGetConverged(eps,&nconv)); | |
| 631 | ✗ | PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of converged eigenpairs: %" PetscInt_FMT "\n\n",nconv)); | |
| 632 | |||
| 633 | ✗ | if (nconv>0) { | |
| 634 | /* Display eigenvalues and relative errors */ | ||
| 635 | ✗ | PetscCall(PetscPrintf(PETSC_COMM_WORLD, | |
| 636 | " k ||Ax-kx||/||kx||\n" | ||
| 637 | " ----------------- ------------------\n")); | ||
| 638 | |||
| 639 | ✗ | for (i=0;i<nconv;i++) { | |
| 640 | /* Get converged eigenpairs: i-th eigenvalue is stored in kr (real part) and | ||
| 641 | ki (imaginary part) */ | ||
| 642 | ✗ | PetscCall(EPSGetEigenpair(eps,i,&kr,&ki,Xr,Xi)); | |
| 643 | /* Compute the relative error associated to each eigenpair */ | ||
| 644 | ✗ | PetscCall(EPSComputeError(eps,i,EPS_ERROR_RELATIVE,&error)); | |
| 645 | |||
| 646 | #if defined(PETSC_USE_COMPLEX) | ||
| 647 | re = PetscRealPart(kr); | ||
| 648 | im = PetscImaginaryPart(kr); | ||
| 649 | #else | ||
| 650 | ✗ | re = kr; | |
| 651 | ✗ | im = ki; | |
| 652 | #endif | ||
| 653 | ✗ | if (im!=0.0) PetscCall(PetscPrintf(PETSC_COMM_WORLD," %9f%+9fi %12g\n",(double)re,(double)im,(double)error)); | |
| 654 | ✗ | else PetscCall(PetscPrintf(PETSC_COMM_WORLD," %12f %12g\n",(double)re,(double)error)); | |
| 655 | } | ||
| 656 | ✗ | PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n")); | |
| 657 | } | ||
| 658 | } | ||
| 659 | |||
| 660 | /* Free work space */ | ||
| 661 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(EPSDestroy(&eps)); |
| 662 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(MatDestroy(&J)); |
| 663 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(MatDestroy(&Jred)); |
| 664 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(MatDestroy(&user.Ybus)); |
| 665 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(VecDestroy(&X)); |
| 666 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(VecDestroy(&Xr)); |
| 667 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(VecDestroy(&Xi)); |
| 668 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(VecDestroy(&user.V0)); |
| 669 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(DMDestroy(&user.dmgen)); |
| 670 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(DMDestroy(&user.dmnet)); |
| 671 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(DMDestroy(&user.dmpgrid)); |
| 672 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(ISDestroy(&user.is_diff)); |
| 673 |
4/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 4 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 1 times.
|
5 | PetscCall(ISDestroy(&user.is_alg)); |
| 674 |
3/6✓ Branch 0 taken 1 times.
✓ Branch 1 taken 1 times.
✓ Branch 2 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
|
5 | PetscCall(SlepcFinalize()); |
| 675 | return 0; | ||
| 676 | } | ||
| 677 | |||
| 678 | /*TEST | ||
| 679 | |||
| 680 | build: | ||
| 681 | requires: !complex | ||
| 682 | |||
| 683 | test: | ||
| 684 | suffix: 1 | ||
| 685 | args: -terse | ||
| 686 | requires: double !complex !defined(PETSC_USE_64BIT_INDICES) | ||
| 687 | localrunfiles: X.bin Ybus.bin | ||
| 688 | |||
| 689 | TEST*/ | ||
| 690 |