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