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 coordiantes.\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 eps = 1.e-7; 326 ierr = VecDuplicate(X,&Y);CHKERRQ(ierr); 327 328 for (i=0;i<ngen;i++) { 329 for (j=0;j<3;j++) PGv[j] = PG[j]; 330 PGv[i] = PG[i]+eps; 331 ierr = InitialGuess(Y,user,PGv);CHKERRQ(ierr); 332 ierr = InitialGuess(X,user,PG);CHKERRQ(ierr); 333 334 ierr = VecAXPY(Y,-1.0,X);CHKERRQ(ierr); 335 ierr = VecScale(Y,1./eps);CHKERRQ(ierr); 336 ierr = VecCopy(Y,DICDP[i]);CHKERRQ(ierr); 337 } 338 ierr = VecDestroy(&Y);CHKERRQ(ierr); 339 PetscFunctionReturn(0); 340 } 341 342 /* Computes F = [-f(x,y);g(x,y)] */ 343 PetscErrorCode ResidualFunction(SNES snes,Vec X, Vec F, Userctx *user) 344 { 345 PetscErrorCode ierr; 346 Vec Xgen,Xnet,Fgen,Fnet; 347 PetscScalar *xgen,*xnet,*fgen,*fnet; 348 PetscInt i,idx=0; 349 PetscScalar Vr,Vi,Vm,Vm2; 350 PetscScalar Eqp,Edp,delta,w; /* Generator variables */ 351 PetscScalar Efd,RF,VR; /* Exciter variables */ 352 PetscScalar Id,Iq; /* Generator dq axis currents */ 353 PetscScalar Vd,Vq,SE; 354 PetscScalar IGr,IGi,IDr,IDi; 355 PetscScalar Zdq_inv[4],det; 356 PetscScalar PD,QD,Vm0,*v0; 357 PetscInt k; 358 359 PetscFunctionBegin; 360 ierr = VecZeroEntries(F);CHKERRQ(ierr); 361 ierr = DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);CHKERRQ(ierr); 362 ierr = DMCompositeGetLocalVectors(user->dmpgrid,&Fgen,&Fnet);CHKERRQ(ierr); 363 ierr = DMCompositeScatter(user->dmpgrid,X,Xgen,Xnet);CHKERRQ(ierr); 364 ierr = DMCompositeScatter(user->dmpgrid,F,Fgen,Fnet);CHKERRQ(ierr); 365 366 /* Network current balance residual IG + Y*V + IL = 0. Only YV is added here. 367 The generator current injection, IG, and load current injection, ID are added later 368 */ 369 /* Note that the values in Ybus are stored assuming the imaginary current balance 370 equation is ordered first followed by real current balance equation for each bus. 371 Thus imaginary current contribution goes in location 2*i, and 372 real current contribution in 2*i+1 373 */ 374 ierr = MatMult(user->Ybus,Xnet,Fnet);CHKERRQ(ierr); 375 376 ierr = VecGetArray(Xgen,&xgen);CHKERRQ(ierr); 377 ierr = VecGetArray(Xnet,&xnet);CHKERRQ(ierr); 378 ierr = VecGetArray(Fgen,&fgen);CHKERRQ(ierr); 379 ierr = VecGetArray(Fnet,&fnet);CHKERRQ(ierr); 380 381 /* Generator subsystem */ 382 for (i=0; i < ngen; i++) { 383 Eqp = xgen[idx]; 384 Edp = xgen[idx+1]; 385 delta = xgen[idx+2]; 386 w = xgen[idx+3]; 387 Id = xgen[idx+4]; 388 Iq = xgen[idx+5]; 389 Efd = xgen[idx+6]; 390 RF = xgen[idx+7]; 391 VR = xgen[idx+8]; 392 393 /* Generator differential equations */ 394 fgen[idx] = (Eqp + (Xd[i] - Xdp[i])*Id - Efd)/Td0p[i]; 395 fgen[idx+1] = (Edp - (Xq[i] - Xqp[i])*Iq)/Tq0p[i]; 396 fgen[idx+2] = -w + w_s; 397 fgen[idx+3] = (-TM[i] + Edp*Id + Eqp*Iq + (Xqp[i] - Xdp[i])*Id*Iq + D[i]*(w - w_s))/M[i]; 398 399 Vr = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */ 400 Vi = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */ 401 402 ierr = ri2dq(Vr,Vi,delta,&Vd,&Vq);CHKERRQ(ierr); 403 /* Algebraic equations for stator currents */ 404 det = Rs[i]*Rs[i] + Xdp[i]*Xqp[i]; 405 406 Zdq_inv[0] = Rs[i]/det; 407 Zdq_inv[1] = Xqp[i]/det; 408 Zdq_inv[2] = -Xdp[i]/det; 409 Zdq_inv[3] = Rs[i]/det; 410 411 fgen[idx+4] = Zdq_inv[0]*(-Edp + Vd) + Zdq_inv[1]*(-Eqp + Vq) + Id; 412 fgen[idx+5] = Zdq_inv[2]*(-Edp + Vd) + Zdq_inv[3]*(-Eqp + Vq) + Iq; 413 414 /* Add generator current injection to network */ 415 ierr = dq2ri(Id,Iq,delta,&IGr,&IGi);CHKERRQ(ierr); 416 417 fnet[2*gbus[i]] -= IGi; 418 fnet[2*gbus[i]+1] -= IGr; 419 420 Vm = PetscSqrtScalar(Vd*Vd + Vq*Vq); 421 422 SE = k1[i]*PetscExpScalar(k2[i]*Efd); 423 424 /* Exciter differential equations */ 425 fgen[idx+6] = (KE[i]*Efd + SE - VR)/TE[i]; 426 fgen[idx+7] = (RF - KF[i]*Efd/TF[i])/TF[i]; 427 fgen[idx+8] = (VR - KA[i]*RF + KA[i]*KF[i]*Efd/TF[i] - KA[i]*(Vref[i] - Vm))/TA[i]; 428 429 idx = idx + 9; 430 } 431 432 ierr = VecGetArray(user->V0,&v0);CHKERRQ(ierr); 433 for (i=0; i < nload; i++) { 434 Vr = xnet[2*lbus[i]]; /* Real part of load bus voltage */ 435 Vi = xnet[2*lbus[i]+1]; /* Imaginary part of the load bus voltage */ 436 Vm = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm; 437 Vm0 = PetscSqrtScalar(v0[2*lbus[i]]*v0[2*lbus[i]] + v0[2*lbus[i]+1]*v0[2*lbus[i]+1]); 438 PD = QD = 0.0; 439 for (k=0; k < ld_nsegsp[i]; k++) PD += ld_alphap[k]*PD0[i]*PetscPowScalar((Vm/Vm0),ld_betap[k]); 440 for (k=0; k < ld_nsegsq[i]; k++) QD += ld_alphaq[k]*QD0[i]*PetscPowScalar((Vm/Vm0),ld_betaq[k]); 441 442 /* Load currents */ 443 IDr = (PD*Vr + QD*Vi)/Vm2; 444 IDi = (-QD*Vr + PD*Vi)/Vm2; 445 446 fnet[2*lbus[i]] += IDi; 447 fnet[2*lbus[i]+1] += IDr; 448 } 449 ierr = VecRestoreArray(user->V0,&v0);CHKERRQ(ierr); 450 451 ierr = VecRestoreArray(Xgen,&xgen);CHKERRQ(ierr); 452 ierr = VecRestoreArray(Xnet,&xnet);CHKERRQ(ierr); 453 ierr = VecRestoreArray(Fgen,&fgen);CHKERRQ(ierr); 454 ierr = VecRestoreArray(Fnet,&fnet);CHKERRQ(ierr); 455 456 ierr = DMCompositeGather(user->dmpgrid,INSERT_VALUES,F,Fgen,Fnet);CHKERRQ(ierr); 457 ierr = DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);CHKERRQ(ierr); 458 ierr = DMCompositeRestoreLocalVectors(user->dmpgrid,&Fgen,&Fnet);CHKERRQ(ierr); 459 PetscFunctionReturn(0); 460 } 461 462 /* \dot{x} - f(x,y) 463 g(x,y) = 0 464 */ 465 PetscErrorCode IFunction(TS ts,PetscReal t, Vec X, Vec Xdot, Vec F, Userctx *user) 466 { 467 PetscErrorCode ierr; 468 SNES snes; 469 PetscScalar *f; 470 const PetscScalar *xdot; 471 PetscInt i; 472 473 PetscFunctionBegin; 474 user->t = t; 475 476 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 477 ierr = ResidualFunction(snes,X,F,user);CHKERRQ(ierr); 478 ierr = VecGetArray(F,&f);CHKERRQ(ierr); 479 ierr = VecGetArrayRead(Xdot,&xdot);CHKERRQ(ierr); 480 for (i=0;i < ngen;i++) { 481 f[9*i] += xdot[9*i]; 482 f[9*i+1] += xdot[9*i+1]; 483 f[9*i+2] += xdot[9*i+2]; 484 f[9*i+3] += xdot[9*i+3]; 485 f[9*i+6] += xdot[9*i+6]; 486 f[9*i+7] += xdot[9*i+7]; 487 f[9*i+8] += xdot[9*i+8]; 488 } 489 ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 490 ierr = VecRestoreArrayRead(Xdot,&xdot);CHKERRQ(ierr); 491 PetscFunctionReturn(0); 492 } 493 494 /* This function is used for solving the algebraic system only during fault on and 495 off times. It computes the entire F and then zeros out the part corresponding to 496 differential equations 497 F = [0;g(y)]; 498 */ 499 PetscErrorCode AlgFunction(SNES snes, Vec X, Vec F, void *ctx) 500 { 501 PetscErrorCode ierr; 502 Userctx *user=(Userctx*)ctx; 503 PetscScalar *f; 504 PetscInt i; 505 506 PetscFunctionBegin; 507 ierr = ResidualFunction(snes,X,F,user);CHKERRQ(ierr); 508 ierr = VecGetArray(F,&f);CHKERRQ(ierr); 509 for (i=0; i < ngen; i++) { 510 f[9*i] = 0; 511 f[9*i+1] = 0; 512 f[9*i+2] = 0; 513 f[9*i+3] = 0; 514 f[9*i+6] = 0; 515 f[9*i+7] = 0; 516 f[9*i+8] = 0; 517 } 518 ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 519 PetscFunctionReturn(0); 520 } 521 522 PetscErrorCode PreallocateJacobian(Mat J, Userctx *user) 523 { 524 PetscErrorCode ierr; 525 PetscInt *d_nnz; 526 PetscInt i,idx=0,start=0; 527 PetscInt ncols; 528 529 PetscFunctionBegin; 530 ierr = PetscMalloc1(user->neqs_pgrid,&d_nnz);CHKERRQ(ierr); 531 for (i=0; i<user->neqs_pgrid; i++) d_nnz[i] = 0; 532 /* Generator subsystem */ 533 for (i=0; i < ngen; i++) { 534 535 d_nnz[idx] += 3; 536 d_nnz[idx+1] += 2; 537 d_nnz[idx+2] += 2; 538 d_nnz[idx+3] += 5; 539 d_nnz[idx+4] += 6; 540 d_nnz[idx+5] += 6; 541 542 d_nnz[user->neqs_gen+2*gbus[i]] += 3; 543 d_nnz[user->neqs_gen+2*gbus[i]+1] += 3; 544 545 d_nnz[idx+6] += 2; 546 d_nnz[idx+7] += 2; 547 d_nnz[idx+8] += 5; 548 549 idx = idx + 9; 550 } 551 552 start = user->neqs_gen; 553 for (i=0; i < nbus; i++) { 554 ierr = MatGetRow(user->Ybus,2*i,&ncols,NULL,NULL);CHKERRQ(ierr); 555 d_nnz[start+2*i] += ncols; 556 d_nnz[start+2*i+1] += ncols; 557 ierr = MatRestoreRow(user->Ybus,2*i,&ncols,NULL,NULL);CHKERRQ(ierr); 558 } 559 560 ierr = MatSeqAIJSetPreallocation(J,0,d_nnz);CHKERRQ(ierr); 561 ierr = PetscFree(d_nnz);CHKERRQ(ierr); 562 PetscFunctionReturn(0); 563 } 564 565 /* 566 J = [-df_dx, -df_dy 567 dg_dx, dg_dy] 568 */ 569 PetscErrorCode ResidualJacobian(SNES snes,Vec X,Mat J,Mat B,void *ctx) 570 { 571 PetscErrorCode ierr; 572 Userctx *user=(Userctx*)ctx; 573 Vec Xgen,Xnet; 574 PetscScalar *xgen,*xnet; 575 PetscInt i,idx=0; 576 PetscScalar Vr,Vi,Vm,Vm2; 577 PetscScalar Eqp,Edp,delta; /* Generator variables */ 578 PetscScalar Efd; /* Exciter variables */ 579 PetscScalar Id,Iq; /* Generator dq axis currents */ 580 PetscScalar Vd,Vq; 581 PetscScalar val[10]; 582 PetscInt row[2],col[10]; 583 PetscInt net_start=user->neqs_gen; 584 PetscInt ncols; 585 const PetscInt *cols; 586 const PetscScalar *yvals; 587 PetscInt k; 588 PetscScalar Zdq_inv[4],det; 589 PetscScalar dVd_dVr,dVd_dVi,dVq_dVr,dVq_dVi,dVd_ddelta,dVq_ddelta; 590 PetscScalar dIGr_ddelta,dIGi_ddelta,dIGr_dId,dIGr_dIq,dIGi_dId,dIGi_dIq; 591 PetscScalar dSE_dEfd; 592 PetscScalar dVm_dVd,dVm_dVq,dVm_dVr,dVm_dVi; 593 PetscScalar PD,QD,Vm0,*v0,Vm4; 594 PetscScalar dPD_dVr,dPD_dVi,dQD_dVr,dQD_dVi; 595 PetscScalar dIDr_dVr,dIDr_dVi,dIDi_dVr,dIDi_dVi; 596 597 PetscFunctionBegin; 598 ierr = MatZeroEntries(B);CHKERRQ(ierr); 599 ierr = DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);CHKERRQ(ierr); 600 ierr = DMCompositeScatter(user->dmpgrid,X,Xgen,Xnet);CHKERRQ(ierr); 601 602 ierr = VecGetArray(Xgen,&xgen);CHKERRQ(ierr); 603 ierr = VecGetArray(Xnet,&xnet);CHKERRQ(ierr); 604 605 /* Generator subsystem */ 606 for (i=0; i < ngen; i++) { 607 Eqp = xgen[idx]; 608 Edp = xgen[idx+1]; 609 delta = xgen[idx+2]; 610 Id = xgen[idx+4]; 611 Iq = xgen[idx+5]; 612 Efd = xgen[idx+6]; 613 614 /* fgen[idx] = (Eqp + (Xd[i] - Xdp[i])*Id - Efd)/Td0p[i]; */ 615 row[0] = idx; 616 col[0] = idx; col[1] = idx+4; col[2] = idx+6; 617 val[0] = 1/ Td0p[i]; val[1] = (Xd[i] - Xdp[i])/ Td0p[i]; val[2] = -1/Td0p[i]; 618 619 ierr = MatSetValues(J,1,row,3,col,val,INSERT_VALUES);CHKERRQ(ierr); 620 621 /* fgen[idx+1] = (Edp - (Xq[i] - Xqp[i])*Iq)/Tq0p[i]; */ 622 row[0] = idx + 1; 623 col[0] = idx + 1; col[1] = idx+5; 624 val[0] = 1/Tq0p[i]; val[1] = -(Xq[i] - Xqp[i])/Tq0p[i]; 625 ierr = MatSetValues(J,1,row,2,col,val,INSERT_VALUES);CHKERRQ(ierr); 626 627 /* fgen[idx+2] = - w + w_s; */ 628 row[0] = idx + 2; 629 col[0] = idx + 2; col[1] = idx + 3; 630 val[0] = 0; val[1] = -1; 631 ierr = MatSetValues(J,1,row,2,col,val,INSERT_VALUES);CHKERRQ(ierr); 632 633 /* fgen[idx+3] = (-TM[i] + Edp*Id + Eqp*Iq + (Xqp[i] - Xdp[i])*Id*Iq + D[i]*(w - w_s))/M[i]; */ 634 row[0] = idx + 3; 635 col[0] = idx; col[1] = idx + 1; col[2] = idx + 3; col[3] = idx + 4; col[4] = idx + 5; 636 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]; 637 ierr = MatSetValues(J,1,row,5,col,val,INSERT_VALUES);CHKERRQ(ierr); 638 639 Vr = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */ 640 Vi = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */ 641 ierr = ri2dq(Vr,Vi,delta,&Vd,&Vq);CHKERRQ(ierr); 642 643 det = Rs[i]*Rs[i] + Xdp[i]*Xqp[i]; 644 645 Zdq_inv[0] = Rs[i]/det; 646 Zdq_inv[1] = Xqp[i]/det; 647 Zdq_inv[2] = -Xdp[i]/det; 648 Zdq_inv[3] = Rs[i]/det; 649 650 dVd_dVr = PetscSinScalar(delta); dVd_dVi = -PetscCosScalar(delta); 651 dVq_dVr = PetscCosScalar(delta); dVq_dVi = PetscSinScalar(delta); 652 dVd_ddelta = Vr*PetscCosScalar(delta) + Vi*PetscSinScalar(delta); 653 dVq_ddelta = -Vr*PetscSinScalar(delta) + Vi*PetscCosScalar(delta); 654 655 /* fgen[idx+4] = Zdq_inv[0]*(-Edp + Vd) + Zdq_inv[1]*(-Eqp + Vq) + Id; */ 656 row[0] = idx+4; 657 col[0] = idx; col[1] = idx+1; col[2] = idx + 2; 658 val[0] = -Zdq_inv[1]; val[1] = -Zdq_inv[0]; val[2] = Zdq_inv[0]*dVd_ddelta + Zdq_inv[1]*dVq_ddelta; 659 col[3] = idx + 4; col[4] = net_start+2*gbus[i]; col[5] = net_start + 2*gbus[i]+1; 660 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; 661 ierr = MatSetValues(J,1,row,6,col,val,INSERT_VALUES);CHKERRQ(ierr); 662 663 /* fgen[idx+5] = Zdq_inv[2]*(-Edp + Vd) + Zdq_inv[3]*(-Eqp + Vq) + Iq; */ 664 row[0] = idx+5; 665 col[0] = idx; col[1] = idx+1; col[2] = idx + 2; 666 val[0] = -Zdq_inv[3]; val[1] = -Zdq_inv[2]; val[2] = Zdq_inv[2]*dVd_ddelta + Zdq_inv[3]*dVq_ddelta; 667 col[3] = idx + 5; col[4] = net_start+2*gbus[i]; col[5] = net_start + 2*gbus[i]+1; 668 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; 669 ierr = MatSetValues(J,1,row,6,col,val,INSERT_VALUES);CHKERRQ(ierr); 670 671 dIGr_ddelta = Id*PetscCosScalar(delta) - Iq*PetscSinScalar(delta); 672 dIGi_ddelta = Id*PetscSinScalar(delta) + Iq*PetscCosScalar(delta); 673 dIGr_dId = PetscSinScalar(delta); dIGr_dIq = PetscCosScalar(delta); 674 dIGi_dId = -PetscCosScalar(delta); dIGi_dIq = PetscSinScalar(delta); 675 676 /* fnet[2*gbus[i]] -= IGi; */ 677 row[0] = net_start + 2*gbus[i]; 678 col[0] = idx+2; col[1] = idx + 4; col[2] = idx + 5; 679 val[0] = -dIGi_ddelta; val[1] = -dIGi_dId; val[2] = -dIGi_dIq; 680 ierr = MatSetValues(J,1,row,3,col,val,INSERT_VALUES);CHKERRQ(ierr); 681 682 /* fnet[2*gbus[i]+1] -= IGr; */ 683 row[0] = net_start + 2*gbus[i]+1; 684 col[0] = idx+2; col[1] = idx + 4; col[2] = idx + 5; 685 val[0] = -dIGr_ddelta; val[1] = -dIGr_dId; val[2] = -dIGr_dIq; 686 ierr = MatSetValues(J,1,row,3,col,val,INSERT_VALUES);CHKERRQ(ierr); 687 688 Vm = PetscSqrtScalar(Vd*Vd + Vq*Vq); 689 690 /* fgen[idx+6] = (KE[i]*Efd + SE - VR)/TE[i]; */ 691 /* SE = k1[i]*PetscExpScalar(k2[i]*Efd); */ 692 dSE_dEfd = k1[i]*k2[i]*PetscExpScalar(k2[i]*Efd); 693 694 row[0] = idx + 6; 695 col[0] = idx + 6; col[1] = idx + 8; 696 val[0] = (KE[i] + dSE_dEfd)/TE[i]; val[1] = -1/TE[i]; 697 ierr = MatSetValues(J,1,row,2,col,val,INSERT_VALUES);CHKERRQ(ierr); 698 699 /* Exciter differential equations */ 700 701 /* fgen[idx+7] = (RF - KF[i]*Efd/TF[i])/TF[i]; */ 702 row[0] = idx + 7; 703 col[0] = idx + 6; col[1] = idx + 7; 704 val[0] = (-KF[i]/TF[i])/TF[i]; val[1] = 1/TF[i]; 705 ierr = MatSetValues(J,1,row,2,col,val,INSERT_VALUES);CHKERRQ(ierr); 706 707 /* fgen[idx+8] = (VR - KA[i]*RF + KA[i]*KF[i]*Efd/TF[i] - KA[i]*(Vref[i] - Vm))/TA[i]; */ 708 /* Vm = (Vd^2 + Vq^2)^0.5; */ 709 dVm_dVd = Vd/Vm; dVm_dVq = Vq/Vm; 710 dVm_dVr = dVm_dVd*dVd_dVr + dVm_dVq*dVq_dVr; 711 dVm_dVi = dVm_dVd*dVd_dVi + dVm_dVq*dVq_dVi; 712 row[0] = idx + 8; 713 col[0] = idx + 6; col[1] = idx + 7; col[2] = idx + 8; 714 val[0] = (KA[i]*KF[i]/TF[i])/TA[i]; val[1] = -KA[i]/TA[i]; val[2] = 1/TA[i]; 715 col[3] = net_start + 2*gbus[i]; col[4] = net_start + 2*gbus[i]+1; 716 val[3] = KA[i]*dVm_dVr/TA[i]; val[4] = KA[i]*dVm_dVi/TA[i]; 717 ierr = MatSetValues(J,1,row,5,col,val,INSERT_VALUES);CHKERRQ(ierr); 718 idx = idx + 9; 719 } 720 721 for (i=0; i<nbus; i++) { 722 ierr = MatGetRow(user->Ybus,2*i,&ncols,&cols,&yvals);CHKERRQ(ierr); 723 row[0] = net_start + 2*i; 724 for (k=0; k<ncols; k++) { 725 col[k] = net_start + cols[k]; 726 val[k] = yvals[k]; 727 } 728 ierr = MatSetValues(J,1,row,ncols,col,val,INSERT_VALUES);CHKERRQ(ierr); 729 ierr = MatRestoreRow(user->Ybus,2*i,&ncols,&cols,&yvals);CHKERRQ(ierr); 730 731 ierr = MatGetRow(user->Ybus,2*i+1,&ncols,&cols,&yvals);CHKERRQ(ierr); 732 row[0] = net_start + 2*i+1; 733 for (k=0; k<ncols; k++) { 734 col[k] = net_start + cols[k]; 735 val[k] = yvals[k]; 736 } 737 ierr = MatSetValues(J,1,row,ncols,col,val,INSERT_VALUES);CHKERRQ(ierr); 738 ierr = MatRestoreRow(user->Ybus,2*i+1,&ncols,&cols,&yvals);CHKERRQ(ierr); 739 } 740 741 ierr = MatAssemblyBegin(J,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr); 742 ierr = MatAssemblyEnd(J,MAT_FLUSH_ASSEMBLY);CHKERRQ(ierr); 743 744 ierr = VecGetArray(user->V0,&v0);CHKERRQ(ierr); 745 for (i=0; i < nload; i++) { 746 Vr = xnet[2*lbus[i]]; /* Real part of load bus voltage */ 747 Vi = xnet[2*lbus[i]+1]; /* Imaginary part of the load bus voltage */ 748 Vm = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2= Vm*Vm; Vm4 = Vm2*Vm2; 749 Vm0 = PetscSqrtScalar(v0[2*lbus[i]]*v0[2*lbus[i]] + v0[2*lbus[i]+1]*v0[2*lbus[i]+1]); 750 PD = QD = 0.0; 751 dPD_dVr = dPD_dVi = dQD_dVr = dQD_dVi = 0.0; 752 for (k=0; k < ld_nsegsp[i]; k++) { 753 PD += ld_alphap[k]*PD0[i]*PetscPowScalar((Vm/Vm0),ld_betap[k]); 754 dPD_dVr += ld_alphap[k]*ld_betap[k]*PD0[i]*PetscPowScalar((1/Vm0),ld_betap[k])*Vr*PetscPowScalar(Vm,(ld_betap[k]-2)); 755 dPD_dVi += ld_alphap[k]*ld_betap[k]*PD0[i]*PetscPowScalar((1/Vm0),ld_betap[k])*Vi*PetscPowScalar(Vm,(ld_betap[k]-2)); 756 } 757 for (k=0; k < ld_nsegsq[i]; k++) { 758 QD += ld_alphaq[k]*QD0[i]*PetscPowScalar((Vm/Vm0),ld_betaq[k]); 759 dQD_dVr += ld_alphaq[k]*ld_betaq[k]*QD0[i]*PetscPowScalar((1/Vm0),ld_betaq[k])*Vr*PetscPowScalar(Vm,(ld_betaq[k]-2)); 760 dQD_dVi += ld_alphaq[k]*ld_betaq[k]*QD0[i]*PetscPowScalar((1/Vm0),ld_betaq[k])*Vi*PetscPowScalar(Vm,(ld_betaq[k]-2)); 761 } 762 763 /* IDr = (PD*Vr + QD*Vi)/Vm2; */ 764 /* IDi = (-QD*Vr + PD*Vi)/Vm2; */ 765 766 dIDr_dVr = (dPD_dVr*Vr + dQD_dVr*Vi + PD)/Vm2 - ((PD*Vr + QD*Vi)*2*Vr)/Vm4; 767 dIDr_dVi = (dPD_dVi*Vr + dQD_dVi*Vi + QD)/Vm2 - ((PD*Vr + QD*Vi)*2*Vi)/Vm4; 768 769 dIDi_dVr = (-dQD_dVr*Vr + dPD_dVr*Vi - QD)/Vm2 - ((-QD*Vr + PD*Vi)*2*Vr)/Vm4; 770 dIDi_dVi = (-dQD_dVi*Vr + dPD_dVi*Vi + PD)/Vm2 - ((-QD*Vr + PD*Vi)*2*Vi)/Vm4; 771 772 /* fnet[2*lbus[i]] += IDi; */ 773 row[0] = net_start + 2*lbus[i]; 774 col[0] = net_start + 2*lbus[i]; col[1] = net_start + 2*lbus[i]+1; 775 val[0] = dIDi_dVr; val[1] = dIDi_dVi; 776 ierr = MatSetValues(J,1,row,2,col,val,ADD_VALUES);CHKERRQ(ierr); 777 /* fnet[2*lbus[i]+1] += IDr; */ 778 row[0] = net_start + 2*lbus[i]+1; 779 col[0] = net_start + 2*lbus[i]; col[1] = net_start + 2*lbus[i]+1; 780 val[0] = dIDr_dVr; val[1] = dIDr_dVi; 781 ierr = MatSetValues(J,1,row,2,col,val,ADD_VALUES);CHKERRQ(ierr); 782 } 783 ierr = VecRestoreArray(user->V0,&v0);CHKERRQ(ierr); 784 785 ierr = VecRestoreArray(Xgen,&xgen);CHKERRQ(ierr); 786 ierr = VecRestoreArray(Xnet,&xnet);CHKERRQ(ierr); 787 788 ierr = DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);CHKERRQ(ierr); 789 790 ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 791 ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 792 PetscFunctionReturn(0); 793 } 794 795 /* 796 J = [I, 0 797 dg_dx, dg_dy] 798 */ 799 PetscErrorCode AlgJacobian(SNES snes,Vec X,Mat A,Mat B,void *ctx) 800 { 801 PetscErrorCode ierr; 802 Userctx *user=(Userctx*)ctx; 803 804 PetscFunctionBegin; 805 ierr = ResidualJacobian(snes,X,A,B,ctx);CHKERRQ(ierr); 806 ierr = MatSetOption(A,MAT_KEEP_NONZERO_PATTERN,PETSC_TRUE);CHKERRQ(ierr); 807 ierr = MatZeroRowsIS(A,user->is_diff,1.0,NULL,NULL);CHKERRQ(ierr); 808 PetscFunctionReturn(0); 809 } 810 811 /* 812 J = [a*I-df_dx, -df_dy 813 dg_dx, dg_dy] 814 */ 815 816 PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,Userctx *user) 817 { 818 PetscErrorCode ierr; 819 SNES snes; 820 PetscScalar atmp = (PetscScalar) a; 821 PetscInt i,row; 822 823 PetscFunctionBegin; 824 user->t = t; 825 826 ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 827 ierr = ResidualJacobian(snes,X,A,B,user);CHKERRQ(ierr); 828 for (i=0;i < ngen;i++) { 829 row = 9*i; 830 ierr = MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);CHKERRQ(ierr); 831 row = 9*i+1; 832 ierr = MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);CHKERRQ(ierr); 833 row = 9*i+2; 834 ierr = MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);CHKERRQ(ierr); 835 row = 9*i+3; 836 ierr = MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);CHKERRQ(ierr); 837 row = 9*i+6; 838 ierr = MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);CHKERRQ(ierr); 839 row = 9*i+7; 840 ierr = MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);CHKERRQ(ierr); 841 row = 9*i+8; 842 ierr = MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);CHKERRQ(ierr); 843 } 844 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 845 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 846 PetscFunctionReturn(0); 847 } 848 849 /* Matrix JacobianP is constant so that it only needs to be evaluated once */ 850 static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A, void *ctx0) 851 { 852 PetscErrorCode ierr; 853 PetscScalar a; 854 PetscInt row,col; 855 Userctx *ctx=(Userctx*)ctx0; 856 857 PetscFunctionBeginUser; 858 859 if (ctx->jacp_flg) { 860 ierr = MatZeroEntries(A);CHKERRQ(ierr); 861 862 for (col=0;col<3;col++) { 863 a = 1.0/M[col]; 864 row = 9*col+3; 865 ierr = MatSetValues(A,1,&row,1,&col,&a,INSERT_VALUES);CHKERRQ(ierr); 866 } 867 868 ctx->jacp_flg = PETSC_FALSE; 869 870 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 871 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 872 } 873 PetscFunctionReturn(0); 874 } 875 876 static PetscErrorCode CostIntegrand(TS ts,PetscReal t,Vec U,Vec R,Userctx *user) 877 { 878 PetscErrorCode ierr; 879 const PetscScalar *u; 880 PetscInt idx; 881 Vec Xgen,Xnet; 882 PetscScalar *r,*xgen; 883 PetscInt i; 884 885 PetscFunctionBegin; 886 ierr = DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);CHKERRQ(ierr); 887 ierr = DMCompositeScatter(user->dmpgrid,U,Xgen,Xnet);CHKERRQ(ierr); 888 889 ierr = VecGetArray(Xgen,&xgen);CHKERRQ(ierr); 890 891 ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 892 ierr = VecGetArray(R,&r);CHKERRQ(ierr); 893 r[0] = 0.; 894 idx = 0; 895 for (i=0;i<ngen;i++) { 896 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); 897 idx += 9; 898 } 899 ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 900 ierr = VecRestoreArray(R,&r);CHKERRQ(ierr); 901 ierr = DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);CHKERRQ(ierr); 902 PetscFunctionReturn(0); 903 } 904 905 static PetscErrorCode DRDUJacobianTranspose(TS ts,PetscReal t,Vec U,Mat DRDU,Mat B,Userctx *user) 906 { 907 PetscErrorCode ierr; 908 Vec Xgen,Xnet,Dgen,Dnet; 909 PetscScalar *xgen,*dgen; 910 PetscInt i; 911 PetscInt idx; 912 Vec drdu_col; 913 PetscScalar *xarr; 914 915 PetscFunctionBegin; 916 ierr = VecDuplicate(U,&drdu_col);CHKERRQ(ierr); 917 ierr = MatDenseGetColumn(DRDU,0,&xarr);CHKERRQ(ierr); 918 ierr = VecPlaceArray(drdu_col,xarr);CHKERRQ(ierr); 919 ierr = DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);CHKERRQ(ierr); 920 ierr = DMCompositeGetLocalVectors(user->dmpgrid,&Dgen,&Dnet);CHKERRQ(ierr); 921 ierr = DMCompositeScatter(user->dmpgrid,U,Xgen,Xnet);CHKERRQ(ierr); 922 ierr = DMCompositeScatter(user->dmpgrid,drdu_col,Dgen,Dnet);CHKERRQ(ierr); 923 924 ierr = VecGetArray(Xgen,&xgen);CHKERRQ(ierr); 925 ierr = VecGetArray(Dgen,&dgen);CHKERRQ(ierr); 926 927 idx = 0; 928 for (i=0;i<ngen;i++) { 929 dgen[idx+3] = 0.; 930 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); 931 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); 932 idx += 9; 933 } 934 935 ierr = VecRestoreArray(Dgen,&dgen);CHKERRQ(ierr); 936 ierr = VecRestoreArray(Xgen,&xgen);CHKERRQ(ierr); 937 ierr = DMCompositeGather(user->dmpgrid,INSERT_VALUES,drdu_col,Dgen,Dnet);CHKERRQ(ierr); 938 ierr = DMCompositeRestoreLocalVectors(user->dmpgrid,&Dgen,&Dnet);CHKERRQ(ierr); 939 ierr = DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);CHKERRQ(ierr); 940 ierr = VecResetArray(drdu_col);CHKERRQ(ierr); 941 ierr = MatDenseRestoreColumn(DRDU,&xarr);CHKERRQ(ierr); 942 ierr = VecDestroy(&drdu_col);CHKERRQ(ierr); 943 PetscFunctionReturn(0); 944 } 945 946 static PetscErrorCode DRDPJacobianTranspose(TS ts,PetscReal t,Vec U,Mat drdp,Userctx *user) 947 { 948 PetscFunctionBegin; 949 PetscFunctionReturn(0); 950 } 951 952 PetscErrorCode ComputeSensiP(Vec lambda,Vec mu,Vec *DICDP,Userctx *user) 953 { 954 PetscErrorCode ierr; 955 PetscScalar *x,*y,sensip; 956 PetscInt i; 957 958 PetscFunctionBegin; 959 ierr = VecGetArray(lambda,&x);CHKERRQ(ierr); 960 ierr = VecGetArray(mu,&y);CHKERRQ(ierr); 961 962 for (i=0;i<3;i++) { 963 ierr = VecDot(lambda,DICDP[i],&sensip);CHKERRQ(ierr); 964 sensip = sensip+y[i]; 965 /* ierr = PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt %D th parameter: %g \n",i,(double)sensip);CHKERRQ(ierr); */ 966 y[i] = sensip; 967 } 968 ierr = VecRestoreArray(mu,&y);CHKERRQ(ierr); 969 PetscFunctionReturn(0); 970 } 971 972 int main(int argc,char **argv) 973 { 974 Userctx user; 975 Vec p; 976 PetscScalar *x_ptr; 977 PetscErrorCode ierr; 978 PetscMPIInt size; 979 PetscInt i; 980 PetscViewer Xview,Ybusview; 981 PetscInt *idx2; 982 Tao tao; 983 KSP ksp; 984 PC pc; 985 Vec lowerb,upperb; 986 987 ierr = PetscInitialize(&argc,&argv,"petscoptions",help);if (ierr) return ierr; 988 ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 989 if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs"); 990 991 user.jacp_flg = PETSC_TRUE; 992 user.neqs_gen = 9*ngen; /* # eqs. for generator subsystem */ 993 user.neqs_net = 2*nbus; /* # eqs. for network subsystem */ 994 user.neqs_pgrid = user.neqs_gen + user.neqs_net; 995 996 /* Create indices for differential and algebraic equations */ 997 ierr = PetscMalloc1(7*ngen,&idx2);CHKERRQ(ierr); 998 for (i=0; i<ngen; i++) { 999 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; 1000 idx2[7*i+4] = 9*i+6; idx2[7*i+5] = 9*i+7; idx2[7*i+6] = 9*i+8; 1001 } 1002 ierr = ISCreateGeneral(PETSC_COMM_WORLD,7*ngen,idx2,PETSC_COPY_VALUES,&user.is_diff);CHKERRQ(ierr); 1003 ierr = ISComplement(user.is_diff,0,user.neqs_pgrid,&user.is_alg);CHKERRQ(ierr); 1004 ierr = PetscFree(idx2);CHKERRQ(ierr); 1005 1006 /* Set run time options */ 1007 ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Transient stability fault options","");CHKERRQ(ierr); 1008 { 1009 user.tfaulton = 1.0; 1010 user.tfaultoff = 1.2; 1011 user.Rfault = 0.0001; 1012 user.faultbus = 8; 1013 ierr = PetscOptionsReal("-tfaulton","","",user.tfaulton,&user.tfaulton,NULL);CHKERRQ(ierr); 1014 ierr = PetscOptionsReal("-tfaultoff","","",user.tfaultoff,&user.tfaultoff,NULL);CHKERRQ(ierr); 1015 ierr = PetscOptionsInt("-faultbus","","",user.faultbus,&user.faultbus,NULL);CHKERRQ(ierr); 1016 user.t0 = 0.0; 1017 user.tmax = 1.3; 1018 ierr = PetscOptionsReal("-t0","","",user.t0,&user.t0,NULL);CHKERRQ(ierr); 1019 ierr = PetscOptionsReal("-tmax","","",user.tmax,&user.tmax,NULL);CHKERRQ(ierr); 1020 user.freq_u = 61.0; 1021 user.freq_l = 59.0; 1022 user.pow = 2; 1023 ierr = PetscOptionsReal("-frequ","","",user.freq_u,&user.freq_u,NULL);CHKERRQ(ierr); 1024 ierr = PetscOptionsReal("-freql","","",user.freq_l,&user.freq_l,NULL);CHKERRQ(ierr); 1025 ierr = PetscOptionsInt("-pow","","",user.pow,&user.pow,NULL);CHKERRQ(ierr); 1026 1027 } 1028 ierr = PetscOptionsEnd();CHKERRQ(ierr); 1029 1030 /* Create DMs for generator and network subsystems */ 1031 ierr = DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,user.neqs_gen,1,1,NULL,&user.dmgen);CHKERRQ(ierr); 1032 ierr = DMSetOptionsPrefix(user.dmgen,"dmgen_");CHKERRQ(ierr); 1033 ierr = DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,user.neqs_net,1,1,NULL,&user.dmnet);CHKERRQ(ierr); 1034 ierr = DMSetOptionsPrefix(user.dmnet,"dmnet_");CHKERRQ(ierr); 1035 ierr = DMSetFromOptions(user.dmnet);CHKERRQ(ierr); 1036 ierr = DMSetUp(user.dmnet);CHKERRQ(ierr); 1037 /* Create a composite DM packer and add the two DMs */ 1038 ierr = DMCompositeCreate(PETSC_COMM_WORLD,&user.dmpgrid);CHKERRQ(ierr); 1039 ierr = DMSetOptionsPrefix(user.dmpgrid,"pgrid_");CHKERRQ(ierr); 1040 ierr = DMSetFromOptions(user.dmgen);CHKERRQ(ierr); 1041 ierr = DMSetUp(user.dmgen);CHKERRQ(ierr); 1042 ierr = DMCompositeAddDM(user.dmpgrid,user.dmgen);CHKERRQ(ierr); 1043 ierr = DMCompositeAddDM(user.dmpgrid,user.dmnet);CHKERRQ(ierr); 1044 1045 /* Read initial voltage vector and Ybus */ 1046 ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,"X.bin",FILE_MODE_READ,&Xview);CHKERRQ(ierr); 1047 ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,"Ybus.bin",FILE_MODE_READ,&Ybusview);CHKERRQ(ierr); 1048 1049 ierr = VecCreate(PETSC_COMM_WORLD,&user.V0);CHKERRQ(ierr); 1050 ierr = VecSetSizes(user.V0,PETSC_DECIDE,user.neqs_net);CHKERRQ(ierr); 1051 ierr = VecLoad(user.V0,Xview);CHKERRQ(ierr); 1052 1053 ierr = MatCreate(PETSC_COMM_WORLD,&user.Ybus);CHKERRQ(ierr); 1054 ierr = MatSetSizes(user.Ybus,PETSC_DECIDE,PETSC_DECIDE,user.neqs_net,user.neqs_net);CHKERRQ(ierr); 1055 ierr = MatSetType(user.Ybus,MATBAIJ);CHKERRQ(ierr); 1056 /* ierr = MatSetBlockSize(ctx->Ybus,2);CHKERRQ(ierr); */ 1057 ierr = MatLoad(user.Ybus,Ybusview);CHKERRQ(ierr); 1058 1059 ierr = PetscViewerDestroy(&Xview);CHKERRQ(ierr); 1060 ierr = PetscViewerDestroy(&Ybusview);CHKERRQ(ierr); 1061 1062 /* Allocate space for Jacobians */ 1063 ierr = MatCreate(PETSC_COMM_WORLD,&user.J);CHKERRQ(ierr); 1064 ierr = MatSetSizes(user.J,PETSC_DECIDE,PETSC_DECIDE,user.neqs_pgrid,user.neqs_pgrid);CHKERRQ(ierr); 1065 ierr = MatSetFromOptions(user.J);CHKERRQ(ierr); 1066 ierr = PreallocateJacobian(user.J,&user);CHKERRQ(ierr); 1067 1068 ierr = MatCreate(PETSC_COMM_WORLD,&user.Jacp);CHKERRQ(ierr); 1069 ierr = MatSetSizes(user.Jacp,PETSC_DECIDE,PETSC_DECIDE,user.neqs_pgrid,3);CHKERRQ(ierr); 1070 ierr = MatSetFromOptions(user.Jacp);CHKERRQ(ierr); 1071 ierr = MatSetUp(user.Jacp);CHKERRQ(ierr); 1072 ierr = MatZeroEntries(user.Jacp);CHKERRQ(ierr); /* initialize to zeros */ 1073 1074 ierr = MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,3,1,NULL,&user.DRDP);CHKERRQ(ierr); 1075 ierr = MatSetUp(user.DRDP);CHKERRQ(ierr); 1076 ierr = MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,user.neqs_pgrid,1,NULL,&user.DRDU);CHKERRQ(ierr); 1077 ierr = MatSetUp(user.DRDU);CHKERRQ(ierr); 1078 1079 /* Create TAO solver and set desired solution method */ 1080 ierr = TaoCreate(PETSC_COMM_WORLD,&tao);CHKERRQ(ierr); 1081 ierr = TaoSetType(tao,TAOBLMVM);CHKERRQ(ierr); 1082 /* 1083 Optimization starts 1084 */ 1085 /* Set initial solution guess */ 1086 ierr = VecCreateSeq(PETSC_COMM_WORLD,3,&p);CHKERRQ(ierr); 1087 ierr = VecGetArray(p,&x_ptr);CHKERRQ(ierr); 1088 x_ptr[0] = PG[0]; x_ptr[1] = PG[1]; x_ptr[2] = PG[2]; 1089 ierr = VecRestoreArray(p,&x_ptr);CHKERRQ(ierr); 1090 1091 ierr = TaoSetInitialVector(tao,p);CHKERRQ(ierr); 1092 /* Set routine for function and gradient evaluation */ 1093 ierr = TaoSetObjectiveAndGradientRoutine(tao,FormFunctionGradient,&user);CHKERRQ(ierr); 1094 1095 /* Set bounds for the optimization */ 1096 ierr = VecDuplicate(p,&lowerb);CHKERRQ(ierr); 1097 ierr = VecDuplicate(p,&upperb);CHKERRQ(ierr); 1098 ierr = VecGetArray(lowerb,&x_ptr);CHKERRQ(ierr); 1099 x_ptr[0] = 0.5; x_ptr[1] = 0.5; x_ptr[2] = 0.5; 1100 ierr = VecRestoreArray(lowerb,&x_ptr);CHKERRQ(ierr); 1101 ierr = VecGetArray(upperb,&x_ptr);CHKERRQ(ierr); 1102 x_ptr[0] = 2.0; x_ptr[1] = 2.0; x_ptr[2] = 2.0; 1103 ierr = VecRestoreArray(upperb,&x_ptr);CHKERRQ(ierr); 1104 ierr = TaoSetVariableBounds(tao,lowerb,upperb);CHKERRQ(ierr); 1105 1106 /* Check for any TAO command line options */ 1107 ierr = TaoSetFromOptions(tao);CHKERRQ(ierr); 1108 ierr = TaoGetKSP(tao,&ksp);CHKERRQ(ierr); 1109 if (ksp) { 1110 ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 1111 ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 1112 } 1113 1114 /* SOLVE THE APPLICATION */ 1115 ierr = TaoSolve(tao);CHKERRQ(ierr); 1116 1117 ierr = VecView(p,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 1118 /* Free TAO data structures */ 1119 ierr = TaoDestroy(&tao);CHKERRQ(ierr); 1120 1121 ierr = DMDestroy(&user.dmgen);CHKERRQ(ierr); 1122 ierr = DMDestroy(&user.dmnet);CHKERRQ(ierr); 1123 ierr = DMDestroy(&user.dmpgrid);CHKERRQ(ierr); 1124 ierr = ISDestroy(&user.is_diff);CHKERRQ(ierr); 1125 ierr = ISDestroy(&user.is_alg);CHKERRQ(ierr); 1126 1127 ierr = MatDestroy(&user.J);CHKERRQ(ierr); 1128 ierr = MatDestroy(&user.Jacp);CHKERRQ(ierr); 1129 ierr = MatDestroy(&user.Ybus);CHKERRQ(ierr); 1130 /* ierr = MatDestroy(&user.Sol);CHKERRQ(ierr); */ 1131 ierr = VecDestroy(&user.V0);CHKERRQ(ierr); 1132 ierr = VecDestroy(&p);CHKERRQ(ierr); 1133 ierr = VecDestroy(&lowerb);CHKERRQ(ierr); 1134 ierr = VecDestroy(&upperb);CHKERRQ(ierr); 1135 ierr = MatDestroy(&user.DRDU);CHKERRQ(ierr); 1136 ierr = MatDestroy(&user.DRDP);CHKERRQ(ierr); 1137 ierr = PetscFinalize(); 1138 return ierr; 1139 } 1140 1141 /* ------------------------------------------------------------------ */ 1142 /* 1143 FormFunction - Evaluates the function and corresponding gradient. 1144 1145 Input Parameters: 1146 tao - the Tao context 1147 X - the input vector 1148 ptr - optional user-defined context, as set by TaoSetObjectiveAndGradientRoutine() 1149 1150 Output Parameters: 1151 f - the newly evaluated function 1152 G - the newly evaluated gradient 1153 */ 1154 PetscErrorCode FormFunctionGradient(Tao tao,Vec P,PetscReal *f,Vec G,void *ctx0) 1155 { 1156 TS ts,quadts; 1157 SNES snes_alg; 1158 PetscErrorCode ierr; 1159 Userctx *ctx = (Userctx*)ctx0; 1160 Vec X; 1161 PetscInt i; 1162 /* sensitivity context */ 1163 PetscScalar *x_ptr; 1164 Vec lambda[1],q; 1165 Vec mu[1]; 1166 PetscInt steps1,steps2,steps3; 1167 Vec DICDP[3]; 1168 Vec F_alg; 1169 PetscInt row_loc,col_loc; 1170 PetscScalar val; 1171 Vec Xdot; 1172 1173 PetscFunctionBegin; 1174 ierr = VecGetArrayRead(P,(const PetscScalar**)&x_ptr);CHKERRQ(ierr); 1175 PG[0] = x_ptr[0]; 1176 PG[1] = x_ptr[1]; 1177 PG[2] = x_ptr[2]; 1178 ierr = VecRestoreArrayRead(P,(const PetscScalar**)&x_ptr);CHKERRQ(ierr); 1179 1180 ctx->stepnum = 0; 1181 1182 ierr = DMCreateGlobalVector(ctx->dmpgrid,&X);CHKERRQ(ierr); 1183 1184 /* Create matrix to save solutions at each time step */ 1185 /* ierr = MatCreateSeqDense(PETSC_COMM_SELF,ctx->neqs_pgrid+1,1002,NULL,&ctx->Sol);CHKERRQ(ierr); */ 1186 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 1187 Create timestepping solver context 1188 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1189 ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr); 1190 ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr); 1191 ierr = TSSetType(ts,TSCN);CHKERRQ(ierr); 1192 ierr = TSSetIFunction(ts,NULL,(TSIFunction) IFunction,ctx);CHKERRQ(ierr); 1193 ierr = TSSetIJacobian(ts,ctx->J,ctx->J,(TSIJacobian)IJacobian,ctx);CHKERRQ(ierr); 1194 ierr = TSSetApplicationContext(ts,ctx);CHKERRQ(ierr); 1195 /* Set RHS JacobianP */ 1196 ierr = TSSetRHSJacobianP(ts,ctx->Jacp,RHSJacobianP,ctx);CHKERRQ(ierr); 1197 1198 ierr = TSCreateQuadratureTS(ts,PETSC_FALSE,&quadts);CHKERRQ(ierr); 1199 ierr = TSSetRHSFunction(quadts,NULL,(TSRHSFunction)CostIntegrand,ctx);CHKERRQ(ierr); 1200 ierr = TSSetRHSJacobian(quadts,ctx->DRDU,ctx->DRDU,(TSRHSJacobian)DRDUJacobianTranspose,ctx);CHKERRQ(ierr); 1201 ierr = TSSetRHSJacobianP(quadts,ctx->DRDP,(TSRHSJacobianP)DRDPJacobianTranspose,ctx);CHKERRQ(ierr); 1202 1203 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 1204 Set initial conditions 1205 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1206 ierr = SetInitialGuess(X,ctx);CHKERRQ(ierr); 1207 1208 /* Approximate DICDP with finite difference, we want to zero out network variables */ 1209 for (i=0;i<3;i++) { 1210 ierr = VecDuplicate(X,&DICDP[i]);CHKERRQ(ierr); 1211 } 1212 ierr = DICDPFiniteDifference(X,DICDP,ctx);CHKERRQ(ierr); 1213 1214 ierr = VecDuplicate(X,&F_alg);CHKERRQ(ierr); 1215 ierr = SNESCreate(PETSC_COMM_WORLD,&snes_alg);CHKERRQ(ierr); 1216 ierr = SNESSetFunction(snes_alg,F_alg,AlgFunction,ctx);CHKERRQ(ierr); 1217 ierr = MatZeroEntries(ctx->J);CHKERRQ(ierr); 1218 ierr = SNESSetJacobian(snes_alg,ctx->J,ctx->J,AlgJacobian,ctx);CHKERRQ(ierr); 1219 ierr = SNESSetOptionsPrefix(snes_alg,"alg_");CHKERRQ(ierr); 1220 ierr = SNESSetFromOptions(snes_alg);CHKERRQ(ierr); 1221 ctx->alg_flg = PETSC_TRUE; 1222 /* Solve the algebraic equations */ 1223 ierr = SNESSolve(snes_alg,NULL,X);CHKERRQ(ierr); 1224 1225 /* Just to set up the Jacobian structure */ 1226 ierr = VecDuplicate(X,&Xdot);CHKERRQ(ierr); 1227 ierr = IJacobian(ts,0.0,X,Xdot,0.0,ctx->J,ctx->J,ctx);CHKERRQ(ierr); 1228 ierr = VecDestroy(&Xdot);CHKERRQ(ierr); 1229 1230 ctx->stepnum++; 1231 1232 /* 1233 Save trajectory of solution so that TSAdjointSolve() may be used 1234 */ 1235 ierr = TSSetSaveTrajectory(ts);CHKERRQ(ierr); 1236 1237 ierr = TSSetTimeStep(ts,0.01);CHKERRQ(ierr); 1238 ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr); 1239 ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 1240 /* ierr = TSSetPostStep(ts,SaveSolution);CHKERRQ(ierr); */ 1241 1242 /* Prefault period */ 1243 ctx->alg_flg = PETSC_FALSE; 1244 ierr = TSSetTime(ts,0.0);CHKERRQ(ierr); 1245 ierr = TSSetMaxTime(ts,ctx->tfaulton);CHKERRQ(ierr); 1246 ierr = TSSolve(ts,X);CHKERRQ(ierr); 1247 ierr = TSGetStepNumber(ts,&steps1);CHKERRQ(ierr); 1248 1249 /* Create the nonlinear solver for solving the algebraic system */ 1250 /* Note that although the algebraic system needs to be solved only for 1251 Idq and V, we reuse the entire system including xgen. The xgen 1252 variables are held constant by setting their residuals to 0 and 1253 putting a 1 on the Jacobian diagonal for xgen rows 1254 */ 1255 ierr = MatZeroEntries(ctx->J);CHKERRQ(ierr); 1256 1257 /* Apply disturbance - resistive fault at ctx->faultbus */ 1258 /* This is done by adding shunt conductance to the diagonal location 1259 in the Ybus matrix */ 1260 row_loc = 2*ctx->faultbus; col_loc = 2*ctx->faultbus+1; /* Location for G */ 1261 val = 1/ctx->Rfault; 1262 ierr = MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);CHKERRQ(ierr); 1263 row_loc = 2*ctx->faultbus+1; col_loc = 2*ctx->faultbus; /* Location for G */ 1264 val = 1/ctx->Rfault; 1265 ierr = MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);CHKERRQ(ierr); 1266 1267 ierr = MatAssemblyBegin(ctx->Ybus,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1268 ierr = MatAssemblyEnd(ctx->Ybus,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1269 1270 ctx->alg_flg = PETSC_TRUE; 1271 /* Solve the algebraic equations */ 1272 ierr = SNESSolve(snes_alg,NULL,X);CHKERRQ(ierr); 1273 1274 ctx->stepnum++; 1275 1276 /* Disturbance period */ 1277 ctx->alg_flg = PETSC_FALSE; 1278 ierr = TSSetTime(ts,ctx->tfaulton);CHKERRQ(ierr); 1279 ierr = TSSetMaxTime(ts,ctx->tfaultoff);CHKERRQ(ierr); 1280 ierr = TSSolve(ts,X);CHKERRQ(ierr); 1281 ierr = TSGetStepNumber(ts,&steps2);CHKERRQ(ierr); 1282 steps2 -= steps1; 1283 1284 /* Remove the fault */ 1285 row_loc = 2*ctx->faultbus; col_loc = 2*ctx->faultbus+1; 1286 val = -1/ctx->Rfault; 1287 ierr = MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);CHKERRQ(ierr); 1288 row_loc = 2*ctx->faultbus+1; col_loc = 2*ctx->faultbus; 1289 val = -1/ctx->Rfault; 1290 ierr = MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);CHKERRQ(ierr); 1291 1292 ierr = MatAssemblyBegin(ctx->Ybus,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1293 ierr = MatAssemblyEnd(ctx->Ybus,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1294 1295 ierr = MatZeroEntries(ctx->J);CHKERRQ(ierr); 1296 1297 ctx->alg_flg = PETSC_TRUE; 1298 1299 /* Solve the algebraic equations */ 1300 ierr = SNESSolve(snes_alg,NULL,X);CHKERRQ(ierr); 1301 1302 ctx->stepnum++; 1303 1304 /* Post-disturbance period */ 1305 ctx->alg_flg = PETSC_TRUE; 1306 ierr = TSSetTime(ts,ctx->tfaultoff);CHKERRQ(ierr); 1307 ierr = TSSetMaxTime(ts,ctx->tmax);CHKERRQ(ierr); 1308 ierr = TSSolve(ts,X);CHKERRQ(ierr); 1309 ierr = TSGetStepNumber(ts,&steps3);CHKERRQ(ierr); 1310 steps3 -= steps2; 1311 steps3 -= steps1; 1312 1313 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 1314 Adjoint model starts here 1315 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1316 ierr = TSSetPostStep(ts,NULL);CHKERRQ(ierr); 1317 ierr = MatCreateVecs(ctx->J,&lambda[0],NULL);CHKERRQ(ierr); 1318 /* Set initial conditions for the adjoint integration */ 1319 ierr = VecZeroEntries(lambda[0]);CHKERRQ(ierr); 1320 1321 ierr = MatCreateVecs(ctx->Jacp,&mu[0],NULL);CHKERRQ(ierr); 1322 ierr = VecZeroEntries(mu[0]);CHKERRQ(ierr); 1323 ierr = TSSetCostGradients(ts,1,lambda,mu);CHKERRQ(ierr); 1324 1325 ierr = TSAdjointSetSteps(ts,steps3);CHKERRQ(ierr); 1326 ierr = TSAdjointSolve(ts);CHKERRQ(ierr); 1327 1328 ierr = MatZeroEntries(ctx->J);CHKERRQ(ierr); 1329 /* Applying disturbance - resistive fault at ctx->faultbus */ 1330 /* This is done by deducting shunt conductance to the diagonal location 1331 in the Ybus matrix */ 1332 row_loc = 2*ctx->faultbus; col_loc = 2*ctx->faultbus+1; /* Location for G */ 1333 val = 1./ctx->Rfault; 1334 ierr = MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);CHKERRQ(ierr); 1335 row_loc = 2*ctx->faultbus+1; col_loc = 2*ctx->faultbus; /* Location for G */ 1336 val = 1./ctx->Rfault; 1337 ierr = MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);CHKERRQ(ierr); 1338 1339 ierr = MatAssemblyBegin(ctx->Ybus,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1340 ierr = MatAssemblyEnd(ctx->Ybus,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1341 1342 /* Set number of steps for the adjoint integration */ 1343 ierr = TSAdjointSetSteps(ts,steps2);CHKERRQ(ierr); 1344 ierr = TSAdjointSolve(ts);CHKERRQ(ierr); 1345 1346 ierr = MatZeroEntries(ctx->J);CHKERRQ(ierr); 1347 /* remove the fault */ 1348 row_loc = 2*ctx->faultbus; col_loc = 2*ctx->faultbus+1; /* Location for G */ 1349 val = -1./ctx->Rfault; 1350 ierr = MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);CHKERRQ(ierr); 1351 row_loc = 2*ctx->faultbus+1; col_loc = 2*ctx->faultbus; /* Location for G */ 1352 val = -1./ctx->Rfault; 1353 ierr = MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);CHKERRQ(ierr); 1354 1355 ierr = MatAssemblyBegin(ctx->Ybus,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1356 ierr = MatAssemblyEnd(ctx->Ybus,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1357 1358 /* Set number of steps for the adjoint integration */ 1359 ierr = TSAdjointSetSteps(ts,steps1);CHKERRQ(ierr); 1360 ierr = TSAdjointSolve(ts);CHKERRQ(ierr); 1361 1362 ierr = ComputeSensiP(lambda[0],mu[0],DICDP,ctx);CHKERRQ(ierr); 1363 ierr = VecCopy(mu[0],G);CHKERRQ(ierr); 1364 1365 ierr = TSGetQuadratureTS(ts,NULL,&quadts);CHKERRQ(ierr); 1366 ierr = TSGetSolution(quadts,&q);CHKERRQ(ierr); 1367 ierr = VecGetArray(q,&x_ptr);CHKERRQ(ierr); 1368 *f = x_ptr[0]; 1369 x_ptr[0] = 0; 1370 ierr = VecRestoreArray(q,&x_ptr);CHKERRQ(ierr); 1371 1372 ierr = VecDestroy(&lambda[0]);CHKERRQ(ierr); 1373 ierr = VecDestroy(&mu[0]);CHKERRQ(ierr); 1374 1375 ierr = SNESDestroy(&snes_alg);CHKERRQ(ierr); 1376 ierr = VecDestroy(&F_alg);CHKERRQ(ierr); 1377 ierr = VecDestroy(&X);CHKERRQ(ierr); 1378 ierr = TSDestroy(&ts);CHKERRQ(ierr); 1379 for (i=0;i<3;i++) { 1380 ierr = VecDestroy(&DICDP[i]);CHKERRQ(ierr); 1381 } 1382 PetscFunctionReturn(0); 1383 } 1384 1385 /*TEST 1386 1387 build: 1388 requires: double !complex !define(PETSC_USE_64BIT_INDICES) 1389 1390 test: 1391 args: -viewer_binary_skip_info -tao_monitor -tao_gttol .2 1392 localrunfiles: petscoptions X.bin Ybus.bin 1393 1394 TEST*/ 1395