Actual source code: ex9busopt.c

  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";

  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: */

 15: #include <petsctao.h>
 16: #include <petscts.h>
 17: #include <petscdm.h>
 18: #include <petscdmda.h>
 19: #include <petscdmcomposite.h>
 20: #include <petsctime.h>

 22: PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal*,Vec,void*);

 24: #define freq 60
 25: #define w_s (2*PETSC_PI*freq)

 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 */

 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};*/

 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 */

 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) */

 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

 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};

 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;

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:   *Fr =  Fd*PetscSinScalar(delta) + Fq*PetscCosScalar(delta);
114:   *Fi = -Fd*PetscCosScalar(delta) + Fq*PetscSinScalar(delta);
115:   return 0;
116: }

118: /* Converts from network frame ([phase a real,imag) to machine (dq) reference frame */
119: PetscErrorCode ri2dq(PetscScalar Fr,PetscScalar Fi,PetscScalar delta,PetscScalar *Fd, PetscScalar *Fq)
120: {
121:   *Fd =  Fr*PetscSinScalar(delta) - Fi*PetscCosScalar(delta);
122:   *Fq =  Fr*PetscCosScalar(delta) + Fi*PetscSinScalar(delta);
123:   return 0;
124: }

126: /* Saves the solution at each time to a matrix */
127: PetscErrorCode SaveSolution(TS ts)
128: {
129:   Userctx           *user;
130:   Vec               X;
131:   PetscScalar       *mat;
132:   const PetscScalar *x;
133:   PetscInt          idx;
134:   PetscReal         t;

136:   TSGetApplicationContext(ts,&user);
137:   TSGetTime(ts,&t);
138:   TSGetSolution(ts,&X);
139:   idx      = user->stepnum*(user->neqs_pgrid+1);
140:   MatDenseGetArray(user->Sol,&mat);
141:   VecGetArrayRead(X,&x);
142:   mat[idx] = t;
143:   PetscArraycpy(mat+idx+1,x,user->neqs_pgrid);
144:   MatDenseRestoreArray(user->Sol,&mat);
145:   VecRestoreArrayRead(X,&x);
146:   user->stepnum++;
147:   return 0;
148: }

150: PetscErrorCode SetInitialGuess(Vec X,Userctx *user)
151: {
152:   Vec            Xgen,Xnet;
153:   PetscScalar    *xgen,*xnet;
154:   PetscInt       i,idx=0;
155:   PetscScalar    Vr,Vi,IGr,IGi,Vm,Vm2;
156:   PetscScalar    Eqp,Edp,delta;
157:   PetscScalar    Efd,RF,VR; /* Exciter variables */
158:   PetscScalar    Id,Iq;  /* Generator dq axis currents */
159:   PetscScalar    theta,Vd,Vq,SE;

161:   M[0] = 2*H[0]/w_s; M[1] = 2*H[1]/w_s; M[2] = 2*H[2]/w_s;
162:   D[0] = 0.1*M[0]; D[1] = 0.1*M[1]; D[2] = 0.1*M[2];

164:   DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);

166:   /* Network subsystem initialization */
167:   VecCopy(user->V0,Xnet);

169:   /* Generator subsystem initialization */
170:   VecGetArray(Xgen,&xgen);
171:   VecGetArray(Xnet,&xnet);

173:   for (i=0; i < ngen; i++) {
174:     Vr  = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */
175:     Vi  = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */
176:     Vm  = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm;
177:     IGr = (Vr*PG[i] + Vi*QG[i])/Vm2;
178:     IGi = (Vi*PG[i] - Vr*QG[i])/Vm2;

180:     delta = PetscAtan2Real(Vi+Xq[i]*IGr,Vr-Xq[i]*IGi); /* Machine angle */

182:     theta = PETSC_PI/2.0 - delta;

184:     Id = IGr*PetscCosScalar(theta) - IGi*PetscSinScalar(theta); /* d-axis stator current */
185:     Iq = IGr*PetscSinScalar(theta) + IGi*PetscCosScalar(theta); /* q-axis stator current */

187:     Vd = Vr*PetscCosScalar(theta) - Vi*PetscSinScalar(theta);
188:     Vq = Vr*PetscSinScalar(theta) + Vi*PetscCosScalar(theta);

190:     Edp = Vd + Rs[i]*Id - Xqp[i]*Iq; /* d-axis transient EMF */
191:     Eqp = Vq + Rs[i]*Iq + Xdp[i]*Id; /* q-axis transient EMF */

193:     TM[i] = PG[i];

195:     /* The generator variables are ordered as [Eqp,Edp,delta,w,Id,Iq] */
196:     xgen[idx]   = Eqp;
197:     xgen[idx+1] = Edp;
198:     xgen[idx+2] = delta;
199:     xgen[idx+3] = w_s;

201:     idx = idx + 4;

203:     xgen[idx]   = Id;
204:     xgen[idx+1] = Iq;

206:     idx = idx + 2;

208:     /* Exciter */
209:     Efd = Eqp + (Xd[i] - Xdp[i])*Id;
210:     SE  = k1[i]*PetscExpScalar(k2[i]*Efd);
211:     VR  =  KE[i]*Efd + SE;
212:     RF  =  KF[i]*Efd/TF[i];

214:     xgen[idx]   = Efd;
215:     xgen[idx+1] = RF;
216:     xgen[idx+2] = VR;

218:     Vref[i] = Vm + (VR/KA[i]);

220:     idx = idx + 3;
221:   }

223:   VecRestoreArray(Xgen,&xgen);
224:   VecRestoreArray(Xnet,&xnet);

226:   /* VecView(Xgen,0); */
227:   DMCompositeGather(user->dmpgrid,INSERT_VALUES,X,Xgen,Xnet);
228:   DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);
229:   return 0;
230: }

232: PetscErrorCode InitialGuess(Vec X,Userctx *user, const PetscScalar PGv[])
233: {
234:   Vec            Xgen,Xnet;
235:   PetscScalar    *xgen,*xnet;
236:   PetscInt       i,idx=0;
237:   PetscScalar    Vr,Vi,IGr,IGi,Vm,Vm2;
238:   PetscScalar    Eqp,Edp,delta;
239:   PetscScalar    Efd,RF,VR; /* Exciter variables */
240:   PetscScalar    Id,Iq;  /* Generator dq axis currents */
241:   PetscScalar    theta,Vd,Vq,SE;

243:   M[0] = 2*H[0]/w_s; M[1] = 2*H[1]/w_s; M[2] = 2*H[2]/w_s;
244:   D[0] = 0.1*M[0]; D[1] = 0.1*M[1]; D[2] = 0.1*M[2];

246:   DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);

248:   /* Network subsystem initialization */
249:   VecCopy(user->V0,Xnet);

251:   /* Generator subsystem initialization */
252:   VecGetArray(Xgen,&xgen);
253:   VecGetArray(Xnet,&xnet);

255:   for (i=0; i < ngen; i++) {
256:     Vr  = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */
257:     Vi  = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */
258:     Vm  = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm;
259:     IGr = (Vr*PGv[i] + Vi*QG[i])/Vm2;
260:     IGi = (Vi*PGv[i] - Vr*QG[i])/Vm2;

262:     delta = PetscAtan2Real(Vi+Xq[i]*IGr,Vr-Xq[i]*IGi); /* Machine angle */

264:     theta = PETSC_PI/2.0 - delta;

266:     Id = IGr*PetscCosScalar(theta) - IGi*PetscSinScalar(theta); /* d-axis stator current */
267:     Iq = IGr*PetscSinScalar(theta) + IGi*PetscCosScalar(theta); /* q-axis stator current */

269:     Vd = Vr*PetscCosScalar(theta) - Vi*PetscSinScalar(theta);
270:     Vq = Vr*PetscSinScalar(theta) + Vi*PetscCosScalar(theta);

272:     Edp = Vd + Rs[i]*Id - Xqp[i]*Iq; /* d-axis transient EMF */
273:     Eqp = Vq + Rs[i]*Iq + Xdp[i]*Id; /* q-axis transient EMF */

275:     /* The generator variables are ordered as [Eqp,Edp,delta,w,Id,Iq] */
276:     xgen[idx]   = Eqp;
277:     xgen[idx+1] = Edp;
278:     xgen[idx+2] = delta;
279:     xgen[idx+3] = w_s;

281:     idx = idx + 4;

283:     xgen[idx]   = Id;
284:     xgen[idx+1] = Iq;

286:     idx = idx + 2;

288:     /* Exciter */
289:     Efd = Eqp + (Xd[i] - Xdp[i])*Id;
290:     SE  = k1[i]*PetscExpScalar(k2[i]*Efd);
291:     VR  =  KE[i]*Efd + SE;
292:     RF  =  KF[i]*Efd/TF[i];

294:     xgen[idx]   = Efd;
295:     xgen[idx+1] = RF;
296:     xgen[idx+2] = VR;

298:     idx = idx + 3;
299:   }

301:   VecRestoreArray(Xgen,&xgen);
302:   VecRestoreArray(Xnet,&xnet);

304:   /* VecView(Xgen,0); */
305:   DMCompositeGather(user->dmpgrid,INSERT_VALUES,X,Xgen,Xnet);
306:   DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);
307:   return 0;
308: }

310: PetscErrorCode DICDPFiniteDifference(Vec X,Vec *DICDP, Userctx *user)
311: {
312:   Vec            Y;
313:   PetscScalar    PGv[3],eps;
314:   PetscInt       i,j;

316:   eps = 1.e-7;
317:   VecDuplicate(X,&Y);

319:   for (i=0;i<ngen;i++) {
320:     for (j=0;j<3;j++) PGv[j] = PG[j];
321:     PGv[i] = PG[i]+eps;
322:     InitialGuess(Y,user,PGv);
323:     InitialGuess(X,user,PG);

325:     VecAXPY(Y,-1.0,X);
326:     VecScale(Y,1./eps);
327:     VecCopy(Y,DICDP[i]);
328:   }
329:   VecDestroy(&Y);
330:   return 0;
331: }

333: /* Computes F = [-f(x,y);g(x,y)] */
334: PetscErrorCode ResidualFunction(SNES snes,Vec X, Vec F, Userctx *user)
335: {
336:   Vec            Xgen,Xnet,Fgen,Fnet;
337:   PetscScalar    *xgen,*xnet,*fgen,*fnet;
338:   PetscInt       i,idx=0;
339:   PetscScalar    Vr,Vi,Vm,Vm2;
340:   PetscScalar    Eqp,Edp,delta,w; /* Generator variables */
341:   PetscScalar    Efd,RF,VR; /* Exciter variables */
342:   PetscScalar    Id,Iq;  /* Generator dq axis currents */
343:   PetscScalar    Vd,Vq,SE;
344:   PetscScalar    IGr,IGi,IDr,IDi;
345:   PetscScalar    Zdq_inv[4],det;
346:   PetscScalar    PD,QD,Vm0,*v0;
347:   PetscInt       k;

349:   VecZeroEntries(F);
350:   DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
351:   DMCompositeGetLocalVectors(user->dmpgrid,&Fgen,&Fnet);
352:   DMCompositeScatter(user->dmpgrid,X,Xgen,Xnet);
353:   DMCompositeScatter(user->dmpgrid,F,Fgen,Fnet);

355:   /* Network current balance residual IG + Y*V + IL = 0. Only YV is added here.
356:      The generator current injection, IG, and load current injection, ID are added later
357:   */
358:   /* Note that the values in Ybus are stored assuming the imaginary current balance
359:      equation is ordered first followed by real current balance equation for each bus.
360:      Thus imaginary current contribution goes in location 2*i, and
361:      real current contribution in 2*i+1
362:   */
363:   MatMult(user->Ybus,Xnet,Fnet);

365:   VecGetArray(Xgen,&xgen);
366:   VecGetArray(Xnet,&xnet);
367:   VecGetArray(Fgen,&fgen);
368:   VecGetArray(Fnet,&fnet);

370:   /* Generator subsystem */
371:   for (i=0; i < ngen; i++) {
372:     Eqp   = xgen[idx];
373:     Edp   = xgen[idx+1];
374:     delta = xgen[idx+2];
375:     w     = xgen[idx+3];
376:     Id    = xgen[idx+4];
377:     Iq    = xgen[idx+5];
378:     Efd   = xgen[idx+6];
379:     RF    = xgen[idx+7];
380:     VR    = xgen[idx+8];

382:     /* Generator differential equations */
383:     fgen[idx]   = (Eqp + (Xd[i] - Xdp[i])*Id - Efd)/Td0p[i];
384:     fgen[idx+1] = (Edp - (Xq[i] - Xqp[i])*Iq)/Tq0p[i];
385:     fgen[idx+2] = -w + w_s;
386:     fgen[idx+3] = (-TM[i] + Edp*Id + Eqp*Iq + (Xqp[i] - Xdp[i])*Id*Iq + D[i]*(w - w_s))/M[i];

388:     Vr = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */
389:     Vi = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */

391:     ri2dq(Vr,Vi,delta,&Vd,&Vq);
392:     /* Algebraic equations for stator currents */
393:     det = Rs[i]*Rs[i] + Xdp[i]*Xqp[i];

395:     Zdq_inv[0] = Rs[i]/det;
396:     Zdq_inv[1] = Xqp[i]/det;
397:     Zdq_inv[2] = -Xdp[i]/det;
398:     Zdq_inv[3] = Rs[i]/det;

400:     fgen[idx+4] = Zdq_inv[0]*(-Edp + Vd) + Zdq_inv[1]*(-Eqp + Vq) + Id;
401:     fgen[idx+5] = Zdq_inv[2]*(-Edp + Vd) + Zdq_inv[3]*(-Eqp + Vq) + Iq;

403:     /* Add generator current injection to network */
404:     dq2ri(Id,Iq,delta,&IGr,&IGi);

406:     fnet[2*gbus[i]]   -= IGi;
407:     fnet[2*gbus[i]+1] -= IGr;

409:     Vm = PetscSqrtScalar(Vd*Vd + Vq*Vq);

411:     SE = k1[i]*PetscExpScalar(k2[i]*Efd);

413:     /* Exciter differential equations */
414:     fgen[idx+6] = (KE[i]*Efd + SE - VR)/TE[i];
415:     fgen[idx+7] = (RF - KF[i]*Efd/TF[i])/TF[i];
416:     fgen[idx+8] = (VR - KA[i]*RF + KA[i]*KF[i]*Efd/TF[i] - KA[i]*(Vref[i] - Vm))/TA[i];

418:     idx = idx + 9;
419:   }

421:   VecGetArray(user->V0,&v0);
422:   for (i=0; i < nload; i++) {
423:     Vr  = xnet[2*lbus[i]]; /* Real part of load bus voltage */
424:     Vi  = xnet[2*lbus[i]+1]; /* Imaginary part of the load bus voltage */
425:     Vm  = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm;
426:     Vm0 = PetscSqrtScalar(v0[2*lbus[i]]*v0[2*lbus[i]] + v0[2*lbus[i]+1]*v0[2*lbus[i]+1]);
427:     PD  = QD = 0.0;
428:     for (k=0; k < ld_nsegsp[i]; k++) PD += ld_alphap[k]*PD0[i]*PetscPowScalar((Vm/Vm0),ld_betap[k]);
429:     for (k=0; k < ld_nsegsq[i]; k++) QD += ld_alphaq[k]*QD0[i]*PetscPowScalar((Vm/Vm0),ld_betaq[k]);

431:     /* Load currents */
432:     IDr = (PD*Vr + QD*Vi)/Vm2;
433:     IDi = (-QD*Vr + PD*Vi)/Vm2;

435:     fnet[2*lbus[i]]   += IDi;
436:     fnet[2*lbus[i]+1] += IDr;
437:   }
438:   VecRestoreArray(user->V0,&v0);

440:   VecRestoreArray(Xgen,&xgen);
441:   VecRestoreArray(Xnet,&xnet);
442:   VecRestoreArray(Fgen,&fgen);
443:   VecRestoreArray(Fnet,&fnet);

445:   DMCompositeGather(user->dmpgrid,INSERT_VALUES,F,Fgen,Fnet);
446:   DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);
447:   DMCompositeRestoreLocalVectors(user->dmpgrid,&Fgen,&Fnet);
448:   return 0;
449: }

451: /* \dot{x} - f(x,y)
452:      g(x,y) = 0
453:  */
454: PetscErrorCode IFunction(TS ts,PetscReal t, Vec X, Vec Xdot, Vec F, Userctx *user)
455: {
456:   SNES              snes;
457:   PetscScalar       *f;
458:   const PetscScalar *xdot;
459:   PetscInt          i;

461:   user->t = t;

463:   TSGetSNES(ts,&snes);
464:   ResidualFunction(snes,X,F,user);
465:   VecGetArray(F,&f);
466:   VecGetArrayRead(Xdot,&xdot);
467:   for (i=0;i < ngen;i++) {
468:     f[9*i]   += xdot[9*i];
469:     f[9*i+1] += xdot[9*i+1];
470:     f[9*i+2] += xdot[9*i+2];
471:     f[9*i+3] += xdot[9*i+3];
472:     f[9*i+6] += xdot[9*i+6];
473:     f[9*i+7] += xdot[9*i+7];
474:     f[9*i+8] += xdot[9*i+8];
475:   }
476:   VecRestoreArray(F,&f);
477:   VecRestoreArrayRead(Xdot,&xdot);
478:   return 0;
479: }

481: /* This function is used for solving the algebraic system only during fault on and
482:    off times. It computes the entire F and then zeros out the part corresponding to
483:    differential equations
484:  F = [0;g(y)];
485: */
486: PetscErrorCode AlgFunction(SNES snes, Vec X, Vec F, void *ctx)
487: {
488:   Userctx        *user=(Userctx*)ctx;
489:   PetscScalar    *f;
490:   PetscInt       i;

492:   ResidualFunction(snes,X,F,user);
493:   VecGetArray(F,&f);
494:   for (i=0; i < ngen; i++) {
495:     f[9*i]   = 0;
496:     f[9*i+1] = 0;
497:     f[9*i+2] = 0;
498:     f[9*i+3] = 0;
499:     f[9*i+6] = 0;
500:     f[9*i+7] = 0;
501:     f[9*i+8] = 0;
502:   }
503:   VecRestoreArray(F,&f);
504:   return 0;
505: }

507: PetscErrorCode PreallocateJacobian(Mat J, Userctx *user)
508: {
509:   PetscInt       *d_nnz;
510:   PetscInt       i,idx=0,start=0;
511:   PetscInt       ncols;

513:   PetscMalloc1(user->neqs_pgrid,&d_nnz);
514:   for (i=0; i<user->neqs_pgrid; i++) d_nnz[i] = 0;
515:   /* Generator subsystem */
516:   for (i=0; i < ngen; i++) {

518:     d_nnz[idx]   += 3;
519:     d_nnz[idx+1] += 2;
520:     d_nnz[idx+2] += 2;
521:     d_nnz[idx+3] += 5;
522:     d_nnz[idx+4] += 6;
523:     d_nnz[idx+5] += 6;

525:     d_nnz[user->neqs_gen+2*gbus[i]]   += 3;
526:     d_nnz[user->neqs_gen+2*gbus[i]+1] += 3;

528:     d_nnz[idx+6] += 2;
529:     d_nnz[idx+7] += 2;
530:     d_nnz[idx+8] += 5;

532:     idx = idx + 9;
533:   }

535:   start = user->neqs_gen;
536:   for (i=0; i < nbus; i++) {
537:     MatGetRow(user->Ybus,2*i,&ncols,NULL,NULL);
538:     d_nnz[start+2*i]   += ncols;
539:     d_nnz[start+2*i+1] += ncols;
540:     MatRestoreRow(user->Ybus,2*i,&ncols,NULL,NULL);
541:   }

543:   MatSeqAIJSetPreallocation(J,0,d_nnz);
544:   PetscFree(d_nnz);
545:   return 0;
546: }

548: /*
549:    J = [-df_dx, -df_dy
550:         dg_dx, dg_dy]
551: */
552: PetscErrorCode ResidualJacobian(SNES snes,Vec X,Mat J,Mat B,void *ctx)
553: {
554:   Userctx           *user=(Userctx*)ctx;
555:   Vec               Xgen,Xnet;
556:   PetscScalar       *xgen,*xnet;
557:   PetscInt          i,idx=0;
558:   PetscScalar       Vr,Vi,Vm,Vm2;
559:   PetscScalar       Eqp,Edp,delta; /* Generator variables */
560:   PetscScalar       Efd; /* Exciter variables */
561:   PetscScalar       Id,Iq;  /* Generator dq axis currents */
562:   PetscScalar       Vd,Vq;
563:   PetscScalar       val[10];
564:   PetscInt          row[2],col[10];
565:   PetscInt          net_start=user->neqs_gen;
566:   PetscInt          ncols;
567:   const PetscInt    *cols;
568:   const PetscScalar *yvals;
569:   PetscInt          k;
570:   PetscScalar       Zdq_inv[4],det;
571:   PetscScalar       dVd_dVr,dVd_dVi,dVq_dVr,dVq_dVi,dVd_ddelta,dVq_ddelta;
572:   PetscScalar       dIGr_ddelta,dIGi_ddelta,dIGr_dId,dIGr_dIq,dIGi_dId,dIGi_dIq;
573:   PetscScalar       dSE_dEfd;
574:   PetscScalar       dVm_dVd,dVm_dVq,dVm_dVr,dVm_dVi;
575:   PetscScalar       PD,QD,Vm0,*v0,Vm4;
576:   PetscScalar       dPD_dVr,dPD_dVi,dQD_dVr,dQD_dVi;
577:   PetscScalar       dIDr_dVr,dIDr_dVi,dIDi_dVr,dIDi_dVi;

579:   MatZeroEntries(B);
580:   DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
581:   DMCompositeScatter(user->dmpgrid,X,Xgen,Xnet);

583:   VecGetArray(Xgen,&xgen);
584:   VecGetArray(Xnet,&xnet);

586:   /* Generator subsystem */
587:   for (i=0; i < ngen; i++) {
588:     Eqp   = xgen[idx];
589:     Edp   = xgen[idx+1];
590:     delta = xgen[idx+2];
591:     Id    = xgen[idx+4];
592:     Iq    = xgen[idx+5];
593:     Efd   = xgen[idx+6];

595:     /*    fgen[idx]   = (Eqp + (Xd[i] - Xdp[i])*Id - Efd)/Td0p[i]; */
596:     row[0] = idx;
597:     col[0] = idx;           col[1] = idx+4;          col[2] = idx+6;
598:     val[0] = 1/ Td0p[i]; val[1] = (Xd[i] - Xdp[i])/ Td0p[i]; val[2] = -1/Td0p[i];

600:     MatSetValues(J,1,row,3,col,val,INSERT_VALUES);

602:     /*    fgen[idx+1] = (Edp - (Xq[i] - Xqp[i])*Iq)/Tq0p[i]; */
603:     row[0] = idx + 1;
604:     col[0] = idx + 1;       col[1] = idx+5;
605:     val[0] = 1/Tq0p[i]; val[1] = -(Xq[i] - Xqp[i])/Tq0p[i];
606:     MatSetValues(J,1,row,2,col,val,INSERT_VALUES);

608:     /*    fgen[idx+2] = - w + w_s; */
609:     row[0] = idx + 2;
610:     col[0] = idx + 2; col[1] = idx + 3;
611:     val[0] = 0;       val[1] = -1;
612:     MatSetValues(J,1,row,2,col,val,INSERT_VALUES);

614:     /*    fgen[idx+3] = (-TM[i] + Edp*Id + Eqp*Iq + (Xqp[i] - Xdp[i])*Id*Iq + D[i]*(w - w_s))/M[i]; */
615:     row[0] = idx + 3;
616:     col[0] = idx; col[1] = idx + 1; col[2] = idx + 3;       col[3] = idx + 4;                  col[4] = idx + 5;
617:     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];
618:     MatSetValues(J,1,row,5,col,val,INSERT_VALUES);

620:     Vr   = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */
621:     Vi   = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */
622:     ri2dq(Vr,Vi,delta,&Vd,&Vq);

624:     det = Rs[i]*Rs[i] + Xdp[i]*Xqp[i];

626:     Zdq_inv[0] = Rs[i]/det;
627:     Zdq_inv[1] = Xqp[i]/det;
628:     Zdq_inv[2] = -Xdp[i]/det;
629:     Zdq_inv[3] = Rs[i]/det;

631:     dVd_dVr    = PetscSinScalar(delta); dVd_dVi = -PetscCosScalar(delta);
632:     dVq_dVr    = PetscCosScalar(delta); dVq_dVi = PetscSinScalar(delta);
633:     dVd_ddelta = Vr*PetscCosScalar(delta) + Vi*PetscSinScalar(delta);
634:     dVq_ddelta = -Vr*PetscSinScalar(delta) + Vi*PetscCosScalar(delta);

636:     /*    fgen[idx+4] = Zdq_inv[0]*(-Edp + Vd) + Zdq_inv[1]*(-Eqp + Vq) + Id; */
637:     row[0] = idx+4;
638:     col[0] = idx;         col[1] = idx+1;        col[2] = idx + 2;
639:     val[0] = -Zdq_inv[1]; val[1] = -Zdq_inv[0];  val[2] = Zdq_inv[0]*dVd_ddelta + Zdq_inv[1]*dVq_ddelta;
640:     col[3] = idx + 4; col[4] = net_start+2*gbus[i];                     col[5] = net_start + 2*gbus[i]+1;
641:     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;
642:     MatSetValues(J,1,row,6,col,val,INSERT_VALUES);

644:     /*  fgen[idx+5] = Zdq_inv[2]*(-Edp + Vd) + Zdq_inv[3]*(-Eqp + Vq) + Iq; */
645:     row[0] = idx+5;
646:     col[0] = idx;         col[1] = idx+1;        col[2] = idx + 2;
647:     val[0] = -Zdq_inv[3]; val[1] = -Zdq_inv[2];  val[2] = Zdq_inv[2]*dVd_ddelta + Zdq_inv[3]*dVq_ddelta;
648:     col[3] = idx + 5; col[4] = net_start+2*gbus[i];                     col[5] = net_start + 2*gbus[i]+1;
649:     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;
650:     MatSetValues(J,1,row,6,col,val,INSERT_VALUES);

652:     dIGr_ddelta = Id*PetscCosScalar(delta) - Iq*PetscSinScalar(delta);
653:     dIGi_ddelta = Id*PetscSinScalar(delta) + Iq*PetscCosScalar(delta);
654:     dIGr_dId    = PetscSinScalar(delta);  dIGr_dIq = PetscCosScalar(delta);
655:     dIGi_dId    = -PetscCosScalar(delta); dIGi_dIq = PetscSinScalar(delta);

657:     /* fnet[2*gbus[i]]   -= IGi; */
658:     row[0] = net_start + 2*gbus[i];
659:     col[0] = idx+2;        col[1] = idx + 4;   col[2] = idx + 5;
660:     val[0] = -dIGi_ddelta; val[1] = -dIGi_dId; val[2] = -dIGi_dIq;
661:     MatSetValues(J,1,row,3,col,val,INSERT_VALUES);

663:     /* fnet[2*gbus[i]+1]   -= IGr; */
664:     row[0] = net_start + 2*gbus[i]+1;
665:     col[0] = idx+2;        col[1] = idx + 4;   col[2] = idx + 5;
666:     val[0] = -dIGr_ddelta; val[1] = -dIGr_dId; val[2] = -dIGr_dIq;
667:     MatSetValues(J,1,row,3,col,val,INSERT_VALUES);

669:     Vm = PetscSqrtScalar(Vd*Vd + Vq*Vq);

671:     /*    fgen[idx+6] = (KE[i]*Efd + SE - VR)/TE[i]; */
672:     /*    SE  = k1[i]*PetscExpScalar(k2[i]*Efd); */
673:     dSE_dEfd = k1[i]*k2[i]*PetscExpScalar(k2[i]*Efd);

675:     row[0] = idx + 6;
676:     col[0] = idx + 6;                     col[1] = idx + 8;
677:     val[0] = (KE[i] + dSE_dEfd)/TE[i];  val[1] = -1/TE[i];
678:     MatSetValues(J,1,row,2,col,val,INSERT_VALUES);

680:     /* Exciter differential equations */

682:     /*    fgen[idx+7] = (RF - KF[i]*Efd/TF[i])/TF[i]; */
683:     row[0] = idx + 7;
684:     col[0] = idx + 6;       col[1] = idx + 7;
685:     val[0] = (-KF[i]/TF[i])/TF[i];  val[1] = 1/TF[i];
686:     MatSetValues(J,1,row,2,col,val,INSERT_VALUES);

688:     /*    fgen[idx+8] = (VR - KA[i]*RF + KA[i]*KF[i]*Efd/TF[i] - KA[i]*(Vref[i] - Vm))/TA[i]; */
689:     /* Vm = (Vd^2 + Vq^2)^0.5; */
690:     dVm_dVd    = Vd/Vm; dVm_dVq = Vq/Vm;
691:     dVm_dVr    = dVm_dVd*dVd_dVr + dVm_dVq*dVq_dVr;
692:     dVm_dVi    = dVm_dVd*dVd_dVi + dVm_dVq*dVq_dVi;
693:     row[0]     = idx + 8;
694:     col[0]     = idx + 6;           col[1] = idx + 7; col[2] = idx + 8;
695:     val[0]     = (KA[i]*KF[i]/TF[i])/TA[i]; val[1] = -KA[i]/TA[i];  val[2] = 1/TA[i];
696:     col[3]     = net_start + 2*gbus[i]; col[4] = net_start + 2*gbus[i]+1;
697:     val[3]     = KA[i]*dVm_dVr/TA[i];         val[4] = KA[i]*dVm_dVi/TA[i];
698:     MatSetValues(J,1,row,5,col,val,INSERT_VALUES);
699:     idx        = idx + 9;
700:   }

702:   for (i=0; i<nbus; i++) {
703:     MatGetRow(user->Ybus,2*i,&ncols,&cols,&yvals);
704:     row[0] = net_start + 2*i;
705:     for (k=0; k<ncols; k++) {
706:       col[k] = net_start + cols[k];
707:       val[k] = yvals[k];
708:     }
709:     MatSetValues(J,1,row,ncols,col,val,INSERT_VALUES);
710:     MatRestoreRow(user->Ybus,2*i,&ncols,&cols,&yvals);

712:     MatGetRow(user->Ybus,2*i+1,&ncols,&cols,&yvals);
713:     row[0] = net_start + 2*i+1;
714:     for (k=0; k<ncols; k++) {
715:       col[k] = net_start + cols[k];
716:       val[k] = yvals[k];
717:     }
718:     MatSetValues(J,1,row,ncols,col,val,INSERT_VALUES);
719:     MatRestoreRow(user->Ybus,2*i+1,&ncols,&cols,&yvals);
720:   }

722:   MatAssemblyBegin(J,MAT_FLUSH_ASSEMBLY);
723:   MatAssemblyEnd(J,MAT_FLUSH_ASSEMBLY);

725:   VecGetArray(user->V0,&v0);
726:   for (i=0; i < nload; i++) {
727:     Vr      = xnet[2*lbus[i]]; /* Real part of load bus voltage */
728:     Vi      = xnet[2*lbus[i]+1]; /* Imaginary part of the load bus voltage */
729:     Vm      = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2= Vm*Vm; Vm4 = Vm2*Vm2;
730:     Vm0     = PetscSqrtScalar(v0[2*lbus[i]]*v0[2*lbus[i]] + v0[2*lbus[i]+1]*v0[2*lbus[i]+1]);
731:     PD      = QD = 0.0;
732:     dPD_dVr = dPD_dVi = dQD_dVr = dQD_dVi = 0.0;
733:     for (k=0; k < ld_nsegsp[i]; k++) {
734:       PD      += ld_alphap[k]*PD0[i]*PetscPowScalar((Vm/Vm0),ld_betap[k]);
735:       dPD_dVr += ld_alphap[k]*ld_betap[k]*PD0[i]*PetscPowScalar((1/Vm0),ld_betap[k])*Vr*PetscPowScalar(Vm,(ld_betap[k]-2));
736:       dPD_dVi += ld_alphap[k]*ld_betap[k]*PD0[i]*PetscPowScalar((1/Vm0),ld_betap[k])*Vi*PetscPowScalar(Vm,(ld_betap[k]-2));
737:     }
738:     for (k=0; k < ld_nsegsq[i]; k++) {
739:       QD      += ld_alphaq[k]*QD0[i]*PetscPowScalar((Vm/Vm0),ld_betaq[k]);
740:       dQD_dVr += ld_alphaq[k]*ld_betaq[k]*QD0[i]*PetscPowScalar((1/Vm0),ld_betaq[k])*Vr*PetscPowScalar(Vm,(ld_betaq[k]-2));
741:       dQD_dVi += ld_alphaq[k]*ld_betaq[k]*QD0[i]*PetscPowScalar((1/Vm0),ld_betaq[k])*Vi*PetscPowScalar(Vm,(ld_betaq[k]-2));
742:     }

744:     /*    IDr = (PD*Vr + QD*Vi)/Vm2; */
745:     /*    IDi = (-QD*Vr + PD*Vi)/Vm2; */

747:     dIDr_dVr = (dPD_dVr*Vr + dQD_dVr*Vi + PD)/Vm2 - ((PD*Vr + QD*Vi)*2*Vr)/Vm4;
748:     dIDr_dVi = (dPD_dVi*Vr + dQD_dVi*Vi + QD)/Vm2 - ((PD*Vr + QD*Vi)*2*Vi)/Vm4;

750:     dIDi_dVr = (-dQD_dVr*Vr + dPD_dVr*Vi - QD)/Vm2 - ((-QD*Vr + PD*Vi)*2*Vr)/Vm4;
751:     dIDi_dVi = (-dQD_dVi*Vr + dPD_dVi*Vi + PD)/Vm2 - ((-QD*Vr + PD*Vi)*2*Vi)/Vm4;

753:     /*    fnet[2*lbus[i]]   += IDi; */
754:     row[0] = net_start + 2*lbus[i];
755:     col[0] = net_start + 2*lbus[i];  col[1] = net_start + 2*lbus[i]+1;
756:     val[0] = dIDi_dVr;               val[1] = dIDi_dVi;
757:     MatSetValues(J,1,row,2,col,val,ADD_VALUES);
758:     /*    fnet[2*lbus[i]+1] += IDr; */
759:     row[0] = net_start + 2*lbus[i]+1;
760:     col[0] = net_start + 2*lbus[i];  col[1] = net_start + 2*lbus[i]+1;
761:     val[0] = dIDr_dVr;               val[1] = dIDr_dVi;
762:     MatSetValues(J,1,row,2,col,val,ADD_VALUES);
763:   }
764:   VecRestoreArray(user->V0,&v0);

766:   VecRestoreArray(Xgen,&xgen);
767:   VecRestoreArray(Xnet,&xnet);

769:   DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);

771:   MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
772:   MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
773:   return 0;
774: }

776: /*
777:    J = [I, 0
778:         dg_dx, dg_dy]
779: */
780: PetscErrorCode AlgJacobian(SNES snes,Vec X,Mat A,Mat B,void *ctx)
781: {
782:   Userctx        *user=(Userctx*)ctx;

784:   ResidualJacobian(snes,X,A,B,ctx);
785:   MatSetOption(A,MAT_KEEP_NONZERO_PATTERN,PETSC_TRUE);
786:   MatZeroRowsIS(A,user->is_diff,1.0,NULL,NULL);
787:   return 0;
788: }

790: /*
791:    J = [a*I-df_dx, -df_dy
792:         dg_dx, dg_dy]
793: */

795: PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,Userctx *user)
796: {
797:   SNES           snes;
798:   PetscScalar    atmp = (PetscScalar) a;
799:   PetscInt       i,row;

801:   user->t = t;

803:   TSGetSNES(ts,&snes);
804:   ResidualJacobian(snes,X,A,B,user);
805:   for (i=0;i < ngen;i++) {
806:     row = 9*i;
807:     MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
808:     row  = 9*i+1;
809:     MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
810:     row  = 9*i+2;
811:     MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
812:     row  = 9*i+3;
813:     MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
814:     row  = 9*i+6;
815:     MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
816:     row  = 9*i+7;
817:     MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
818:     row  = 9*i+8;
819:     MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
820:   }
821:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
822:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
823:   return 0;
824: }

826: /* Matrix JacobianP is constant so that it only needs to be evaluated once */
827: static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A, void *ctx0)
828: {
829:   PetscScalar    a;
830:   PetscInt       row,col;
831:   Userctx        *ctx=(Userctx*)ctx0;


835:   if (ctx->jacp_flg) {
836:     MatZeroEntries(A);

838:     for (col=0;col<3;col++) {
839:       a    = 1.0/M[col];
840:       row  = 9*col+3;
841:       MatSetValues(A,1,&row,1,&col,&a,INSERT_VALUES);
842:     }

844:     ctx->jacp_flg = PETSC_FALSE;

846:     MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
847:     MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
848:   }
849:   return 0;
850: }

852: static PetscErrorCode CostIntegrand(TS ts,PetscReal t,Vec U,Vec R,Userctx *user)
853: {
854:   const PetscScalar *u;
855:   PetscInt          idx;
856:   Vec               Xgen,Xnet;
857:   PetscScalar       *r,*xgen;
858:   PetscInt          i;

860:   DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
861:   DMCompositeScatter(user->dmpgrid,U,Xgen,Xnet);

863:   VecGetArray(Xgen,&xgen);

865:   VecGetArrayRead(U,&u);
866:   VecGetArray(R,&r);
867:   r[0] = 0.;
868:   idx = 0;
869:   for (i=0;i<ngen;i++) {
870:     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);
871:     idx  += 9;
872:   }
873:   VecRestoreArrayRead(U,&u);
874:   VecRestoreArray(R,&r);
875:   DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);
876:   return 0;
877: }

879: static PetscErrorCode DRDUJacobianTranspose(TS ts,PetscReal t,Vec U,Mat DRDU,Mat B,Userctx *user)
880: {
881:   Vec            Xgen,Xnet,Dgen,Dnet;
882:   PetscScalar    *xgen,*dgen;
883:   PetscInt       i;
884:   PetscInt       idx;
885:   Vec            drdu_col;
886:   PetscScalar    *xarr;

888:   VecDuplicate(U,&drdu_col);
889:   MatDenseGetColumn(DRDU,0,&xarr);
890:   VecPlaceArray(drdu_col,xarr);
891:   DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
892:   DMCompositeGetLocalVectors(user->dmpgrid,&Dgen,&Dnet);
893:   DMCompositeScatter(user->dmpgrid,U,Xgen,Xnet);
894:   DMCompositeScatter(user->dmpgrid,drdu_col,Dgen,Dnet);

896:   VecGetArray(Xgen,&xgen);
897:   VecGetArray(Dgen,&dgen);

899:   idx = 0;
900:   for (i=0;i<ngen;i++) {
901:     dgen[idx+3] = 0.;
902:     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);
903:     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);
904:     idx += 9;
905:   }

907:   VecRestoreArray(Dgen,&dgen);
908:   VecRestoreArray(Xgen,&xgen);
909:   DMCompositeGather(user->dmpgrid,INSERT_VALUES,drdu_col,Dgen,Dnet);
910:   DMCompositeRestoreLocalVectors(user->dmpgrid,&Dgen,&Dnet);
911:   DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);
912:   VecResetArray(drdu_col);
913:   MatDenseRestoreColumn(DRDU,&xarr);
914:   VecDestroy(&drdu_col);
915:   return 0;
916: }

918: static PetscErrorCode DRDPJacobianTranspose(TS ts,PetscReal t,Vec U,Mat drdp,Userctx *user)
919: {
920:   return 0;
921: }

923: PetscErrorCode ComputeSensiP(Vec lambda,Vec mu,Vec *DICDP,Userctx *user)
924: {
925:   PetscScalar    *x,*y,sensip;
926:   PetscInt       i;

928:   VecGetArray(lambda,&x);
929:   VecGetArray(mu,&y);

931:   for (i=0;i<3;i++) {
932:     VecDot(lambda,DICDP[i],&sensip);
933:     sensip = sensip+y[i];
934:     /* PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt %D th parameter: %g \n",i,(double)sensip); */
935:      y[i] = sensip;
936:   }
937:   VecRestoreArray(mu,&y);
938:   return 0;
939: }

941: int main(int argc,char **argv)
942: {
943:   Userctx            user;
944:   Vec                p;
945:   PetscScalar        *x_ptr;
946:   PetscErrorCode     ierr;
947:   PetscMPIInt        size;
948:   PetscInt           i;
949:   PetscViewer        Xview,Ybusview;
950:   PetscInt           *idx2;
951:   Tao                tao;
952:   KSP                ksp;
953:   PC                 pc;
954:   Vec                lowerb,upperb;

956:   PetscInitialize(&argc,&argv,"petscoptions",help);
957:   MPI_Comm_size(PETSC_COMM_WORLD,&size);

960:   user.jacp_flg   = PETSC_TRUE;
961:   user.neqs_gen   = 9*ngen; /* # eqs. for generator subsystem */
962:   user.neqs_net   = 2*nbus; /* # eqs. for network subsystem   */
963:   user.neqs_pgrid = user.neqs_gen + user.neqs_net;

965:   /* Create indices for differential and algebraic equations */
966:   PetscMalloc1(7*ngen,&idx2);
967:   for (i=0; i<ngen; i++) {
968:     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;
969:     idx2[7*i+4] = 9*i+6; idx2[7*i+5] = 9*i+7; idx2[7*i+6] = 9*i+8;
970:   }
971:   ISCreateGeneral(PETSC_COMM_WORLD,7*ngen,idx2,PETSC_COPY_VALUES,&user.is_diff);
972:   ISComplement(user.is_diff,0,user.neqs_pgrid,&user.is_alg);
973:   PetscFree(idx2);

975:   /* Set run time options */
976:   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Transient stability fault options","");
977:   {
978:     user.tfaulton  = 1.0;
979:     user.tfaultoff = 1.2;
980:     user.Rfault    = 0.0001;
981:     user.faultbus  = 8;
982:     PetscOptionsReal("-tfaulton","","",user.tfaulton,&user.tfaulton,NULL);
983:     PetscOptionsReal("-tfaultoff","","",user.tfaultoff,&user.tfaultoff,NULL);
984:     PetscOptionsInt("-faultbus","","",user.faultbus,&user.faultbus,NULL);
985:     user.t0        = 0.0;
986:     user.tmax      = 1.3;
987:     PetscOptionsReal("-t0","","",user.t0,&user.t0,NULL);
988:     PetscOptionsReal("-tmax","","",user.tmax,&user.tmax,NULL);
989:     user.freq_u    = 61.0;
990:     user.freq_l    = 59.0;
991:     user.pow       = 2;
992:     PetscOptionsReal("-frequ","","",user.freq_u,&user.freq_u,NULL);
993:     PetscOptionsReal("-freql","","",user.freq_l,&user.freq_l,NULL);
994:     PetscOptionsInt("-pow","","",user.pow,&user.pow,NULL);

996:   }
997:   PetscOptionsEnd();

999:   /* Create DMs for generator and network subsystems */
1000:   DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,user.neqs_gen,1,1,NULL,&user.dmgen);
1001:   DMSetOptionsPrefix(user.dmgen,"dmgen_");
1002:   DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,user.neqs_net,1,1,NULL,&user.dmnet);
1003:   DMSetOptionsPrefix(user.dmnet,"dmnet_");
1004:   DMSetFromOptions(user.dmnet);
1005:   DMSetUp(user.dmnet);
1006:   /* Create a composite DM packer and add the two DMs */
1007:   DMCompositeCreate(PETSC_COMM_WORLD,&user.dmpgrid);
1008:   DMSetOptionsPrefix(user.dmpgrid,"pgrid_");
1009:   DMSetFromOptions(user.dmgen);
1010:   DMSetUp(user.dmgen);
1011:   DMCompositeAddDM(user.dmpgrid,user.dmgen);
1012:   DMCompositeAddDM(user.dmpgrid,user.dmnet);

1014:   /* Read initial voltage vector and Ybus */
1015:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"X.bin",FILE_MODE_READ,&Xview);
1016:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"Ybus.bin",FILE_MODE_READ,&Ybusview);

1018:   VecCreate(PETSC_COMM_WORLD,&user.V0);
1019:   VecSetSizes(user.V0,PETSC_DECIDE,user.neqs_net);
1020:   VecLoad(user.V0,Xview);

1022:   MatCreate(PETSC_COMM_WORLD,&user.Ybus);
1023:   MatSetSizes(user.Ybus,PETSC_DECIDE,PETSC_DECIDE,user.neqs_net,user.neqs_net);
1024:   MatSetType(user.Ybus,MATBAIJ);
1025:   /*  MatSetBlockSize(ctx->Ybus,2); */
1026:   MatLoad(user.Ybus,Ybusview);

1028:   PetscViewerDestroy(&Xview);
1029:   PetscViewerDestroy(&Ybusview);

1031:   /* Allocate space for Jacobians */
1032:   MatCreate(PETSC_COMM_WORLD,&user.J);
1033:   MatSetSizes(user.J,PETSC_DECIDE,PETSC_DECIDE,user.neqs_pgrid,user.neqs_pgrid);
1034:   MatSetFromOptions(user.J);
1035:   PreallocateJacobian(user.J,&user);

1037:   MatCreate(PETSC_COMM_WORLD,&user.Jacp);
1038:   MatSetSizes(user.Jacp,PETSC_DECIDE,PETSC_DECIDE,user.neqs_pgrid,3);
1039:   MatSetFromOptions(user.Jacp);
1040:   MatSetUp(user.Jacp);
1041:   MatZeroEntries(user.Jacp); /* initialize to zeros */

1043:   MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,3,1,NULL,&user.DRDP);
1044:   MatSetUp(user.DRDP);
1045:   MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,user.neqs_pgrid,1,NULL,&user.DRDU);
1046:   MatSetUp(user.DRDU);

1048:   /* Create TAO solver and set desired solution method */
1049:   TaoCreate(PETSC_COMM_WORLD,&tao);
1050:   TaoSetType(tao,TAOBLMVM);
1051:   /*
1052:      Optimization starts
1053:   */
1054:   /* Set initial solution guess */
1055:   VecCreateSeq(PETSC_COMM_WORLD,3,&p);
1056:   VecGetArray(p,&x_ptr);
1057:   x_ptr[0] = PG[0]; x_ptr[1] = PG[1]; x_ptr[2] = PG[2];
1058:   VecRestoreArray(p,&x_ptr);

1060:   TaoSetSolution(tao,p);
1061:   /* Set routine for function and gradient evaluation */
1062:   TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,&user);

1064:   /* Set bounds for the optimization */
1065:   VecDuplicate(p,&lowerb);
1066:   VecDuplicate(p,&upperb);
1067:   VecGetArray(lowerb,&x_ptr);
1068:   x_ptr[0] = 0.5; x_ptr[1] = 0.5; x_ptr[2] = 0.5;
1069:   VecRestoreArray(lowerb,&x_ptr);
1070:   VecGetArray(upperb,&x_ptr);
1071:   x_ptr[0] = 2.0; x_ptr[1] = 2.0; x_ptr[2] = 2.0;
1072:   VecRestoreArray(upperb,&x_ptr);
1073:   TaoSetVariableBounds(tao,lowerb,upperb);

1075:   /* Check for any TAO command line options */
1076:   TaoSetFromOptions(tao);
1077:   TaoGetKSP(tao,&ksp);
1078:   if (ksp) {
1079:     KSPGetPC(ksp,&pc);
1080:     PCSetType(pc,PCNONE);
1081:   }

1083:   /* SOLVE THE APPLICATION */
1084:   TaoSolve(tao);

1086:   VecView(p,PETSC_VIEWER_STDOUT_WORLD);
1087:   /* Free TAO data structures */
1088:   TaoDestroy(&tao);

1090:   DMDestroy(&user.dmgen);
1091:   DMDestroy(&user.dmnet);
1092:   DMDestroy(&user.dmpgrid);
1093:   ISDestroy(&user.is_diff);
1094:   ISDestroy(&user.is_alg);

1096:   MatDestroy(&user.J);
1097:   MatDestroy(&user.Jacp);
1098:   MatDestroy(&user.Ybus);
1099:   /* MatDestroy(&user.Sol); */
1100:   VecDestroy(&user.V0);
1101:   VecDestroy(&p);
1102:   VecDestroy(&lowerb);
1103:   VecDestroy(&upperb);
1104:   MatDestroy(&user.DRDU);
1105:   MatDestroy(&user.DRDP);
1106:   PetscFinalize();
1107:   return 0;
1108: }

1110: /* ------------------------------------------------------------------ */
1111: /*
1112:    FormFunction - Evaluates the function and corresponding gradient.

1114:    Input Parameters:
1115:    tao - the Tao context
1116:    X   - the input vector
1117:    ptr - optional user-defined context, as set by TaoSetObjectiveAndGradient()

1119:    Output Parameters:
1120:    f   - the newly evaluated function
1121:    G   - the newly evaluated gradient
1122: */
1123: PetscErrorCode FormFunctionGradient(Tao tao,Vec P,PetscReal *f,Vec G,void *ctx0)
1124: {
1125:   TS             ts,quadts;
1126:   SNES           snes_alg;
1127:   Userctx        *ctx = (Userctx*)ctx0;
1128:   Vec            X;
1129:   PetscInt       i;
1130:   /* sensitivity context */
1131:   PetscScalar    *x_ptr;
1132:   Vec            lambda[1],q;
1133:   Vec            mu[1];
1134:   PetscInt       steps1,steps2,steps3;
1135:   Vec            DICDP[3];
1136:   Vec            F_alg;
1137:   PetscInt       row_loc,col_loc;
1138:   PetscScalar    val;
1139:   Vec            Xdot;

1141:   VecGetArrayRead(P,(const PetscScalar**)&x_ptr);
1142:   PG[0] = x_ptr[0];
1143:   PG[1] = x_ptr[1];
1144:   PG[2] = x_ptr[2];
1145:   VecRestoreArrayRead(P,(const PetscScalar**)&x_ptr);

1147:   ctx->stepnum = 0;

1149:   DMCreateGlobalVector(ctx->dmpgrid,&X);

1151:   /* Create matrix to save solutions at each time step */
1152:   /* MatCreateSeqDense(PETSC_COMM_SELF,ctx->neqs_pgrid+1,1002,NULL,&ctx->Sol); */
1153:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1154:      Create timestepping solver context
1155:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1156:   TSCreate(PETSC_COMM_WORLD,&ts);
1157:   TSSetProblemType(ts,TS_NONLINEAR);
1158:   TSSetType(ts,TSCN);
1159:   TSSetIFunction(ts,NULL,(TSIFunction) IFunction,ctx);
1160:   TSSetIJacobian(ts,ctx->J,ctx->J,(TSIJacobian)IJacobian,ctx);
1161:   TSSetApplicationContext(ts,ctx);
1162:   /*   Set RHS JacobianP */
1163:   TSSetRHSJacobianP(ts,ctx->Jacp,RHSJacobianP,ctx);

1165:   TSCreateQuadratureTS(ts,PETSC_FALSE,&quadts);
1166:   TSSetRHSFunction(quadts,NULL,(TSRHSFunction)CostIntegrand,ctx);
1167:   TSSetRHSJacobian(quadts,ctx->DRDU,ctx->DRDU,(TSRHSJacobian)DRDUJacobianTranspose,ctx);
1168:   TSSetRHSJacobianP(quadts,ctx->DRDP,(TSRHSJacobianP)DRDPJacobianTranspose,ctx);

1170:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1171:      Set initial conditions
1172:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1173:   SetInitialGuess(X,ctx);

1175:   /* Approximate DICDP with finite difference, we want to zero out network variables */
1176:   for (i=0;i<3;i++) {
1177:     VecDuplicate(X,&DICDP[i]);
1178:   }
1179:   DICDPFiniteDifference(X,DICDP,ctx);

1181:   VecDuplicate(X,&F_alg);
1182:   SNESCreate(PETSC_COMM_WORLD,&snes_alg);
1183:   SNESSetFunction(snes_alg,F_alg,AlgFunction,ctx);
1184:   MatZeroEntries(ctx->J);
1185:   SNESSetJacobian(snes_alg,ctx->J,ctx->J,AlgJacobian,ctx);
1186:   SNESSetOptionsPrefix(snes_alg,"alg_");
1187:   SNESSetFromOptions(snes_alg);
1188:   ctx->alg_flg = PETSC_TRUE;
1189:   /* Solve the algebraic equations */
1190:   SNESSolve(snes_alg,NULL,X);

1192:   /* Just to set up the Jacobian structure */
1193:   VecDuplicate(X,&Xdot);
1194:   IJacobian(ts,0.0,X,Xdot,0.0,ctx->J,ctx->J,ctx);
1195:   VecDestroy(&Xdot);

1197:   ctx->stepnum++;

1199:   /*
1200:     Save trajectory of solution so that TSAdjointSolve() may be used
1201:   */
1202:   TSSetSaveTrajectory(ts);

1204:   TSSetTimeStep(ts,0.01);
1205:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);
1206:   TSSetFromOptions(ts);
1207:   /* TSSetPostStep(ts,SaveSolution); */

1209:   /* Prefault period */
1210:   ctx->alg_flg = PETSC_FALSE;
1211:   TSSetTime(ts,0.0);
1212:   TSSetMaxTime(ts,ctx->tfaulton);
1213:   TSSolve(ts,X);
1214:   TSGetStepNumber(ts,&steps1);

1216:   /* Create the nonlinear solver for solving the algebraic system */
1217:   /* Note that although the algebraic system needs to be solved only for
1218:      Idq and V, we reuse the entire system including xgen. The xgen
1219:      variables are held constant by setting their residuals to 0 and
1220:      putting a 1 on the Jacobian diagonal for xgen rows
1221:   */
1222:   MatZeroEntries(ctx->J);

1224:   /* Apply disturbance - resistive fault at ctx->faultbus */
1225:   /* This is done by adding shunt conductance to the diagonal location
1226:      in the Ybus matrix */
1227:   row_loc = 2*ctx->faultbus; col_loc = 2*ctx->faultbus+1; /* Location for G */
1228:   val     = 1/ctx->Rfault;
1229:   MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);
1230:   row_loc = 2*ctx->faultbus+1; col_loc = 2*ctx->faultbus; /* Location for G */
1231:   val     = 1/ctx->Rfault;
1232:   MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);

1234:   MatAssemblyBegin(ctx->Ybus,MAT_FINAL_ASSEMBLY);
1235:   MatAssemblyEnd(ctx->Ybus,MAT_FINAL_ASSEMBLY);

1237:   ctx->alg_flg = PETSC_TRUE;
1238:   /* Solve the algebraic equations */
1239:   SNESSolve(snes_alg,NULL,X);

1241:   ctx->stepnum++;

1243:   /* Disturbance period */
1244:   ctx->alg_flg = PETSC_FALSE;
1245:   TSSetTime(ts,ctx->tfaulton);
1246:   TSSetMaxTime(ts,ctx->tfaultoff);
1247:   TSSolve(ts,X);
1248:   TSGetStepNumber(ts,&steps2);
1249:   steps2 -= steps1;

1251:   /* Remove the fault */
1252:   row_loc = 2*ctx->faultbus; col_loc = 2*ctx->faultbus+1;
1253:   val     = -1/ctx->Rfault;
1254:   MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);
1255:   row_loc = 2*ctx->faultbus+1; col_loc = 2*ctx->faultbus;
1256:   val     = -1/ctx->Rfault;
1257:   MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);

1259:   MatAssemblyBegin(ctx->Ybus,MAT_FINAL_ASSEMBLY);
1260:   MatAssemblyEnd(ctx->Ybus,MAT_FINAL_ASSEMBLY);

1262:   MatZeroEntries(ctx->J);

1264:   ctx->alg_flg = PETSC_TRUE;

1266:   /* Solve the algebraic equations */
1267:   SNESSolve(snes_alg,NULL,X);

1269:   ctx->stepnum++;

1271:   /* Post-disturbance period */
1272:   ctx->alg_flg = PETSC_TRUE;
1273:   TSSetTime(ts,ctx->tfaultoff);
1274:   TSSetMaxTime(ts,ctx->tmax);
1275:   TSSolve(ts,X);
1276:   TSGetStepNumber(ts,&steps3);
1277:   steps3 -= steps2;
1278:   steps3 -= steps1;

1280:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1281:      Adjoint model starts here
1282:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1283:   TSSetPostStep(ts,NULL);
1284:   MatCreateVecs(ctx->J,&lambda[0],NULL);
1285:   /*   Set initial conditions for the adjoint integration */
1286:   VecZeroEntries(lambda[0]);

1288:   MatCreateVecs(ctx->Jacp,&mu[0],NULL);
1289:   VecZeroEntries(mu[0]);
1290:   TSSetCostGradients(ts,1,lambda,mu);

1292:   TSAdjointSetSteps(ts,steps3);
1293:   TSAdjointSolve(ts);

1295:   MatZeroEntries(ctx->J);
1296:   /* Applying disturbance - resistive fault at ctx->faultbus */
1297:   /* This is done by deducting shunt conductance to the diagonal location
1298:      in the Ybus matrix */
1299:   row_loc = 2*ctx->faultbus; col_loc = 2*ctx->faultbus+1; /* Location for G */
1300:   val     = 1./ctx->Rfault;
1301:   MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);
1302:   row_loc = 2*ctx->faultbus+1; col_loc = 2*ctx->faultbus; /* Location for G */
1303:   val     = 1./ctx->Rfault;
1304:   MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);

1306:   MatAssemblyBegin(ctx->Ybus,MAT_FINAL_ASSEMBLY);
1307:   MatAssemblyEnd(ctx->Ybus,MAT_FINAL_ASSEMBLY);

1309:   /*   Set number of steps for the adjoint integration */
1310:   TSAdjointSetSteps(ts,steps2);
1311:   TSAdjointSolve(ts);

1313:   MatZeroEntries(ctx->J);
1314:   /* remove the fault */
1315:   row_loc = 2*ctx->faultbus; col_loc = 2*ctx->faultbus+1; /* Location for G */
1316:   val     = -1./ctx->Rfault;
1317:   MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);
1318:   row_loc = 2*ctx->faultbus+1; col_loc = 2*ctx->faultbus; /* Location for G */
1319:   val     = -1./ctx->Rfault;
1320:   MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);

1322:   MatAssemblyBegin(ctx->Ybus,MAT_FINAL_ASSEMBLY);
1323:   MatAssemblyEnd(ctx->Ybus,MAT_FINAL_ASSEMBLY);

1325:   /*   Set number of steps for the adjoint integration */
1326:   TSAdjointSetSteps(ts,steps1);
1327:   TSAdjointSolve(ts);

1329:   ComputeSensiP(lambda[0],mu[0],DICDP,ctx);
1330:   VecCopy(mu[0],G);

1332:   TSGetQuadratureTS(ts,NULL,&quadts);
1333:   TSGetSolution(quadts,&q);
1334:   VecGetArray(q,&x_ptr);
1335:   *f   = x_ptr[0];
1336:   x_ptr[0] = 0;
1337:   VecRestoreArray(q,&x_ptr);

1339:   VecDestroy(&lambda[0]);
1340:   VecDestroy(&mu[0]);

1342:   SNESDestroy(&snes_alg);
1343:   VecDestroy(&F_alg);
1344:   VecDestroy(&X);
1345:   TSDestroy(&ts);
1346:   for (i=0;i<3;i++) {
1347:     VecDestroy(&DICDP[i]);
1348:   }
1349:   return 0;
1350: }

1352: /*TEST

1354:    build:
1355:       requires: double !complex !defined(PETSC_USE_64BIT_INDICES)

1357:    test:
1358:       args: -viewer_binary_skip_info -tao_monitor -tao_gttol .2
1359:       localrunfiles: petscoptions X.bin Ybus.bin

1361: TEST*/