Actual source code: ex9busopt.c

petsc-3.9.1 2018-04-29
Report Typos and Errors
  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";

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

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: {
123:   *Fd =  Fr*PetscSinScalar(delta) - Fi*PetscCosScalar(delta);
124:   *Fq =  Fr*PetscCosScalar(delta) + Fi*PetscSinScalar(delta);
125:   return(0);
126: }

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;

140:   TSGetApplicationContext(ts,&user);
141:   TSGetTime(ts,&t);
142:   TSGetSolution(ts,&X);
143:   idx      = user->stepnum*(user->neqs_pgrid+1);
144:   MatDenseGetArray(user->Sol,&mat);
145:   VecGetArrayRead(X,&x);
146:   mat[idx] = t;
147:   PetscMemcpy(mat+idx+1,x,user->neqs_pgrid*sizeof(PetscScalar));
148:   MatDenseRestoreArray(user->Sol,&mat);
149:   VecRestoreArrayRead(X,&x);
150:   user->stepnum++;
151:   return(0);
152: }

154: PetscErrorCode SetInitialGuess(Vec X,Userctx *user)
155: {
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;

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

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

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

175:   /* Generator subsystem initialization */
176:   VecGetArray(Xgen,&xgen);
177:   VecGetArray(Xnet,&xnet);

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;

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

188:     theta = PETSC_PI/2.0 - delta;

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

193:     Vd = Vr*PetscCosScalar(theta) - Vi*PetscSinScalar(theta);
194:     Vq = Vr*PetscSinScalar(theta) + Vi*PetscCosScalar(theta);

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

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

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;

207:     idx = idx + 4;

209:     xgen[idx]   = Id;
210:     xgen[idx+1] = Iq;

212:     idx = idx + 2;

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

220:     xgen[idx]   = Efd;
221:     xgen[idx+1] = RF;
222:     xgen[idx+2] = VR;

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

226:     idx = idx + 3;
227:   }

229:   VecRestoreArray(Xgen,&xgen);
230:   VecRestoreArray(Xnet,&xnet);

232:   /* VecView(Xgen,0); */
233:   DMCompositeGather(user->dmpgrid,INSERT_VALUES,X,Xgen,Xnet);
234:   DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);
235:   return(0);
236: }

238: PetscErrorCode InitialGuess(Vec X,Userctx *user, const PetscScalar PGv[])
239: {
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;

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

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

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

259:   /* Generator subsystem initialization */
260:   VecGetArray(Xgen,&xgen);
261:   VecGetArray(Xnet,&xnet);

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;

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

272:     theta = PETSC_PI/2.0 - delta;

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

277:     Vd = Vr*PetscCosScalar(theta) - Vi*PetscSinScalar(theta);
278:     Vq = Vr*PetscSinScalar(theta) + Vi*PetscCosScalar(theta);

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

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;

289:     idx = idx + 4;

291:     xgen[idx]   = Id;
292:     xgen[idx+1] = Iq;

294:     idx = idx + 2;

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

302:     xgen[idx]   = Efd;
303:     xgen[idx+1] = RF;
304:     xgen[idx+2] = VR;

306:     idx = idx + 3;
307:   }

309:   VecRestoreArray(Xgen,&xgen);
310:   VecRestoreArray(Xnet,&xnet);

312:   /* VecView(Xgen,0); */
313:   DMCompositeGather(user->dmpgrid,INSERT_VALUES,X,Xgen,Xnet);
314:   DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);
315:   return(0);
316: }

318: PetscErrorCode DICDPFiniteDifference(Vec X,Vec *DICDP, Userctx *user)
319: {
320:   Vec            Y;
321:   PetscScalar    PGv[3],eps;
323:   PetscInt       i,j;

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

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:     InitialGuess(Y,user,PGv);
332:     InitialGuess(X,user,PG);

334:     VecAXPY(Y,-1.0,X);
335:     VecScale(Y,1./eps);
336:     VecCopy(Y,DICDP[i]);
337:   }
338:   VecDestroy(&Y);
339:   return(0);
340: }


343: /* Computes F = [-f(x,y);g(x,y)] */
344: PetscErrorCode ResidualFunction(SNES snes,Vec X, Vec F, Userctx *user)
345: {
347:   Vec            Xgen,Xnet,Fgen,Fnet;
348:   PetscScalar    *xgen,*xnet,*fgen,*fnet;
349:   PetscInt       i,idx=0;
350:   PetscScalar    Vr,Vi,Vm,Vm2;
351:   PetscScalar    Eqp,Edp,delta,w; /* Generator variables */
352:   PetscScalar    Efd,RF,VR; /* Exciter variables */
353:   PetscScalar    Id,Iq;  /* Generator dq axis currents */
354:   PetscScalar    Vd,Vq,SE;
355:   PetscScalar    IGr,IGi,IDr,IDi;
356:   PetscScalar    Zdq_inv[4],det;
357:   PetscScalar    PD,QD,Vm0,*v0;
358:   PetscInt       k;

361:   VecZeroEntries(F);
362:   DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
363:   DMCompositeGetLocalVectors(user->dmpgrid,&Fgen,&Fnet);
364:   DMCompositeScatter(user->dmpgrid,X,Xgen,Xnet);
365:   DMCompositeScatter(user->dmpgrid,F,Fgen,Fnet);

367:   /* Network current balance residual IG + Y*V + IL = 0. Only YV is added here.
368:      The generator current injection, IG, and load current injection, ID are added later
369:   */
370:   /* Note that the values in Ybus are stored assuming the imaginary current balance
371:      equation is ordered first followed by real current balance equation for each bus.
372:      Thus imaginary current contribution goes in location 2*i, and
373:      real current contribution in 2*i+1
374:   */
375:   MatMult(user->Ybus,Xnet,Fnet);

377:   VecGetArray(Xgen,&xgen);
378:   VecGetArray(Xnet,&xnet);
379:   VecGetArray(Fgen,&fgen);
380:   VecGetArray(Fnet,&fnet);

382:   /* Generator subsystem */
383:   for (i=0; i < ngen; i++) {
384:     Eqp   = xgen[idx];
385:     Edp   = xgen[idx+1];
386:     delta = xgen[idx+2];
387:     w     = xgen[idx+3];
388:     Id    = xgen[idx+4];
389:     Iq    = xgen[idx+5];
390:     Efd   = xgen[idx+6];
391:     RF    = xgen[idx+7];
392:     VR    = xgen[idx+8];

394:     /* Generator differential equations */
395:     fgen[idx]   = (Eqp + (Xd[i] - Xdp[i])*Id - Efd)/Td0p[i];
396:     fgen[idx+1] = (Edp - (Xq[i] - Xqp[i])*Iq)/Tq0p[i];
397:     fgen[idx+2] = -w + w_s;
398:     fgen[idx+3] = (-TM[i] + Edp*Id + Eqp*Iq + (Xqp[i] - Xdp[i])*Id*Iq + D[i]*(w - w_s))/M[i];

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

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

407:     Zdq_inv[0] = Rs[i]/det;
408:     Zdq_inv[1] = Xqp[i]/det;
409:     Zdq_inv[2] = -Xdp[i]/det;
410:     Zdq_inv[3] = Rs[i]/det;

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

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

418:     fnet[2*gbus[i]]   -= IGi;
419:     fnet[2*gbus[i]+1] -= IGr;

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

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

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

430:     idx = idx + 9;
431:   }

433:   VecGetArray(user->V0,&v0);
434:   for (i=0; i < nload; i++) {
435:     Vr  = xnet[2*lbus[i]]; /* Real part of load bus voltage */
436:     Vi  = xnet[2*lbus[i]+1]; /* Imaginary part of the load bus voltage */
437:     Vm  = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm;
438:     Vm0 = PetscSqrtScalar(v0[2*lbus[i]]*v0[2*lbus[i]] + v0[2*lbus[i]+1]*v0[2*lbus[i]+1]);
439:     PD  = QD = 0.0;
440:     for (k=0; k < ld_nsegsp[i]; k++) PD += ld_alphap[k]*PD0[i]*PetscPowScalar((Vm/Vm0),ld_betap[k]);
441:     for (k=0; k < ld_nsegsq[i]; k++) QD += ld_alphaq[k]*QD0[i]*PetscPowScalar((Vm/Vm0),ld_betaq[k]);

443:     /* Load currents */
444:     IDr = (PD*Vr + QD*Vi)/Vm2;
445:     IDi = (-QD*Vr + PD*Vi)/Vm2;

447:     fnet[2*lbus[i]]   += IDi;
448:     fnet[2*lbus[i]+1] += IDr;
449:   }
450:   VecRestoreArray(user->V0,&v0);

452:   VecRestoreArray(Xgen,&xgen);
453:   VecRestoreArray(Xnet,&xnet);
454:   VecRestoreArray(Fgen,&fgen);
455:   VecRestoreArray(Fnet,&fnet);

457:   DMCompositeGather(user->dmpgrid,INSERT_VALUES,F,Fgen,Fnet);
458:   DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);
459:   DMCompositeRestoreLocalVectors(user->dmpgrid,&Fgen,&Fnet);
460:   return(0);
461: }

463: /* \dot{x} - f(x,y)
464:      g(x,y) = 0
465:  */
466: PetscErrorCode IFunction(TS ts,PetscReal t, Vec X, Vec Xdot, Vec F, Userctx *user)
467: {
468:   PetscErrorCode    ierr;
469:   SNES              snes;
470:   PetscScalar       *f;
471:   const PetscScalar *xdot;
472:   PetscInt          i;

475:   user->t = t;

477:   TSGetSNES(ts,&snes);
478:   ResidualFunction(snes,X,F,user);
479:   VecGetArray(F,&f);
480:   VecGetArrayRead(Xdot,&xdot);
481:   for (i=0;i < ngen;i++) {
482:     f[9*i]   += xdot[9*i];
483:     f[9*i+1] += xdot[9*i+1];
484:     f[9*i+2] += xdot[9*i+2];
485:     f[9*i+3] += xdot[9*i+3];
486:     f[9*i+6] += xdot[9*i+6];
487:     f[9*i+7] += xdot[9*i+7];
488:     f[9*i+8] += xdot[9*i+8];
489:   }
490:   VecRestoreArray(F,&f);
491:   VecRestoreArrayRead(Xdot,&xdot);
492:   return(0);
493: }

495: /* This function is used for solving the algebraic system only during fault on and
496:    off times. It computes the entire F and then zeros out the part corresponding to
497:    differential equations
498:  F = [0;g(y)];
499: */
500: PetscErrorCode AlgFunction(SNES snes, Vec X, Vec F, void *ctx)
501: {
503:   Userctx        *user=(Userctx*)ctx;
504:   PetscScalar    *f;
505:   PetscInt       i;

508:   ResidualFunction(snes,X,F,user);
509:   VecGetArray(F,&f);
510:   for (i=0; i < ngen; i++) {
511:     f[9*i]   = 0;
512:     f[9*i+1] = 0;
513:     f[9*i+2] = 0;
514:     f[9*i+3] = 0;
515:     f[9*i+6] = 0;
516:     f[9*i+7] = 0;
517:     f[9*i+8] = 0;
518:   }
519:   VecRestoreArray(F,&f);
520:   return(0);
521: }

523: PetscErrorCode PreallocateJacobian(Mat J, Userctx *user)
524: {
526:   PetscInt       *d_nnz;
527:   PetscInt       i,idx=0,start=0;
528:   PetscInt       ncols;

531:   PetscMalloc1(user->neqs_pgrid,&d_nnz);
532:   for (i=0; i<user->neqs_pgrid; i++) d_nnz[i] = 0;
533:   /* Generator subsystem */
534:   for (i=0; i < ngen; i++) {

536:     d_nnz[idx]   += 3;
537:     d_nnz[idx+1] += 2;
538:     d_nnz[idx+2] += 2;
539:     d_nnz[idx+3] += 5;
540:     d_nnz[idx+4] += 6;
541:     d_nnz[idx+5] += 6;

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

546:     d_nnz[idx+6] += 2;
547:     d_nnz[idx+7] += 2;
548:     d_nnz[idx+8] += 5;

550:     idx = idx + 9;
551:   }

553:   start = user->neqs_gen;
554:   for (i=0; i < nbus; i++) {
555:     MatGetRow(user->Ybus,2*i,&ncols,NULL,NULL);
556:     d_nnz[start+2*i]   += ncols;
557:     d_nnz[start+2*i+1] += ncols;
558:     MatRestoreRow(user->Ybus,2*i,&ncols,NULL,NULL);
559:   }

561:   MatSeqAIJSetPreallocation(J,0,d_nnz);
562:   PetscFree(d_nnz);
563:   return(0);
564: }

566: /*
567:    J = [-df_dx, -df_dy
568:         dg_dx, dg_dy]
569: */
570: PetscErrorCode ResidualJacobian(SNES snes,Vec X,Mat J,Mat B,void *ctx)
571: {
572:   PetscErrorCode    ierr;
573:   Userctx           *user=(Userctx*)ctx;
574:   Vec               Xgen,Xnet;
575:   PetscScalar       *xgen,*xnet;
576:   PetscInt          i,idx=0;
577:   PetscScalar       Vr,Vi,Vm,Vm2;
578:   PetscScalar       Eqp,Edp,delta; /* Generator variables */
579:   PetscScalar       Efd; /* Exciter variables */
580:   PetscScalar       Id,Iq;  /* Generator dq axis currents */
581:   PetscScalar       Vd,Vq;
582:   PetscScalar       val[10];
583:   PetscInt          row[2],col[10];
584:   PetscInt          net_start=user->neqs_gen;
585:   PetscInt          ncols;
586:   const PetscInt    *cols;
587:   const PetscScalar *yvals;
588:   PetscInt          k;
589:   PetscScalar       Zdq_inv[4],det;
590:   PetscScalar       dVd_dVr,dVd_dVi,dVq_dVr,dVq_dVi,dVd_ddelta,dVq_ddelta;
591:   PetscScalar       dIGr_ddelta,dIGi_ddelta,dIGr_dId,dIGr_dIq,dIGi_dId,dIGi_dIq;
592:   PetscScalar       dSE_dEfd;
593:   PetscScalar       dVm_dVd,dVm_dVq,dVm_dVr,dVm_dVi;
594:   PetscScalar       PD,QD,Vm0,*v0,Vm4;
595:   PetscScalar       dPD_dVr,dPD_dVi,dQD_dVr,dQD_dVi;
596:   PetscScalar       dIDr_dVr,dIDr_dVi,dIDi_dVr,dIDi_dVi;

599:   MatZeroEntries(B);
600:   DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
601:   DMCompositeScatter(user->dmpgrid,X,Xgen,Xnet);

603:   VecGetArray(Xgen,&xgen);
604:   VecGetArray(Xnet,&xnet);

606:   /* Generator subsystem */
607:   for (i=0; i < ngen; i++) {
608:     Eqp   = xgen[idx];
609:     Edp   = xgen[idx+1];
610:     delta = xgen[idx+2];
611:     Id    = xgen[idx+4];
612:     Iq    = xgen[idx+5];
613:     Efd   = xgen[idx+6];

615:     /*    fgen[idx]   = (Eqp + (Xd[i] - Xdp[i])*Id - Efd)/Td0p[i]; */
616:     row[0] = idx;
617:     col[0] = idx;           col[1] = idx+4;          col[2] = idx+6;
618:     val[0] = 1/ Td0p[i]; val[1] = (Xd[i] - Xdp[i])/ Td0p[i]; val[2] = -1/Td0p[i];

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

622:     /*    fgen[idx+1] = (Edp - (Xq[i] - Xqp[i])*Iq)/Tq0p[i]; */
623:     row[0] = idx + 1;
624:     col[0] = idx + 1;       col[1] = idx+5;
625:     val[0] = 1/Tq0p[i]; val[1] = -(Xq[i] - Xqp[i])/Tq0p[i];
626:     MatSetValues(J,1,row,2,col,val,INSERT_VALUES);

628:     /*    fgen[idx+2] = - w + w_s; */
629:     row[0] = idx + 2;
630:     col[0] = idx + 2; col[1] = idx + 3;
631:     val[0] = 0;       val[1] = -1;
632:     MatSetValues(J,1,row,2,col,val,INSERT_VALUES);

634:     /*    fgen[idx+3] = (-TM[i] + Edp*Id + Eqp*Iq + (Xqp[i] - Xdp[i])*Id*Iq + D[i]*(w - w_s))/M[i]; */
635:     row[0] = idx + 3;
636:     col[0] = idx; col[1] = idx + 1; col[2] = idx + 3;       col[3] = idx + 4;                  col[4] = idx + 5;
637:     val[0] = Iq/M[i];  val[1] = Id/M[i];      val[2] = D[i]/M[i]; val[3] = (Edp + (Xqp[i]-Xdp[i])*Iq)/M[i]; val[4] = (Eqp + (Xqp[i] - Xdp[i])*Id)/M[i];
638:     MatSetValues(J,1,row,5,col,val,INSERT_VALUES);

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

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

646:     Zdq_inv[0] = Rs[i]/det;
647:     Zdq_inv[1] = Xqp[i]/det;
648:     Zdq_inv[2] = -Xdp[i]/det;
649:     Zdq_inv[3] = Rs[i]/det;

651:     dVd_dVr    = PetscSinScalar(delta); dVd_dVi = -PetscCosScalar(delta);
652:     dVq_dVr    = PetscCosScalar(delta); dVq_dVi = PetscSinScalar(delta);
653:     dVd_ddelta = Vr*PetscCosScalar(delta) + Vi*PetscSinScalar(delta);
654:     dVq_ddelta = -Vr*PetscSinScalar(delta) + Vi*PetscCosScalar(delta);

656:     /*    fgen[idx+4] = Zdq_inv[0]*(-Edp + Vd) + Zdq_inv[1]*(-Eqp + Vq) + Id; */
657:     row[0] = idx+4;
658:     col[0] = idx;         col[1] = idx+1;        col[2] = idx + 2;
659:     val[0] = -Zdq_inv[1]; val[1] = -Zdq_inv[0];  val[2] = Zdq_inv[0]*dVd_ddelta + Zdq_inv[1]*dVq_ddelta;
660:     col[3] = idx + 4; col[4] = net_start+2*gbus[i];                     col[5] = net_start + 2*gbus[i]+1;
661:     val[3] = 1;       val[4] = Zdq_inv[0]*dVd_dVr + Zdq_inv[1]*dVq_dVr; val[5] = Zdq_inv[0]*dVd_dVi + Zdq_inv[1]*dVq_dVi;
662:     MatSetValues(J,1,row,6,col,val,INSERT_VALUES);

664:     /*  fgen[idx+5] = Zdq_inv[2]*(-Edp + Vd) + Zdq_inv[3]*(-Eqp + Vq) + Iq; */
665:     row[0] = idx+5;
666:     col[0] = idx;         col[1] = idx+1;        col[2] = idx + 2;
667:     val[0] = -Zdq_inv[3]; val[1] = -Zdq_inv[2];  val[2] = Zdq_inv[2]*dVd_ddelta + Zdq_inv[3]*dVq_ddelta;
668:     col[3] = idx + 5; col[4] = net_start+2*gbus[i];                     col[5] = net_start + 2*gbus[i]+1;
669:     val[3] = 1;       val[4] = Zdq_inv[2]*dVd_dVr + Zdq_inv[3]*dVq_dVr; val[5] = Zdq_inv[2]*dVd_dVi + Zdq_inv[3]*dVq_dVi;
670:     MatSetValues(J,1,row,6,col,val,INSERT_VALUES);

672:     dIGr_ddelta = Id*PetscCosScalar(delta) - Iq*PetscSinScalar(delta);
673:     dIGi_ddelta = Id*PetscSinScalar(delta) + Iq*PetscCosScalar(delta);
674:     dIGr_dId    = PetscSinScalar(delta);  dIGr_dIq = PetscCosScalar(delta);
675:     dIGi_dId    = -PetscCosScalar(delta); dIGi_dIq = PetscSinScalar(delta);

677:     /* fnet[2*gbus[i]]   -= IGi; */
678:     row[0] = net_start + 2*gbus[i];
679:     col[0] = idx+2;        col[1] = idx + 4;   col[2] = idx + 5;
680:     val[0] = -dIGi_ddelta; val[1] = -dIGi_dId; val[2] = -dIGi_dIq;
681:     MatSetValues(J,1,row,3,col,val,INSERT_VALUES);

683:     /* fnet[2*gbus[i]+1]   -= IGr; */
684:     row[0] = net_start + 2*gbus[i]+1;
685:     col[0] = idx+2;        col[1] = idx + 4;   col[2] = idx + 5;
686:     val[0] = -dIGr_ddelta; val[1] = -dIGr_dId; val[2] = -dIGr_dIq;
687:     MatSetValues(J,1,row,3,col,val,INSERT_VALUES);

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

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

695:     row[0] = idx + 6;
696:     col[0] = idx + 6;                     col[1] = idx + 8;
697:     val[0] = (KE[i] + dSE_dEfd)/TE[i];  val[1] = -1/TE[i];
698:     MatSetValues(J,1,row,2,col,val,INSERT_VALUES);

700:     /* Exciter differential equations */

702:     /*    fgen[idx+7] = (RF - KF[i]*Efd/TF[i])/TF[i]; */
703:     row[0] = idx + 7;
704:     col[0] = idx + 6;       col[1] = idx + 7;
705:     val[0] = (-KF[i]/TF[i])/TF[i];  val[1] = 1/TF[i];
706:     MatSetValues(J,1,row,2,col,val,INSERT_VALUES);

708:     /*    fgen[idx+8] = (VR - KA[i]*RF + KA[i]*KF[i]*Efd/TF[i] - KA[i]*(Vref[i] - Vm))/TA[i]; */
709:     /* Vm = (Vd^2 + Vq^2)^0.5; */
710:     dVm_dVd    = Vd/Vm; dVm_dVq = Vq/Vm;
711:     dVm_dVr    = dVm_dVd*dVd_dVr + dVm_dVq*dVq_dVr;
712:     dVm_dVi    = dVm_dVd*dVd_dVi + dVm_dVq*dVq_dVi;
713:     row[0]     = idx + 8;
714:     col[0]     = idx + 6;           col[1] = idx + 7; col[2] = idx + 8;
715:     val[0]     = (KA[i]*KF[i]/TF[i])/TA[i]; val[1] = -KA[i]/TA[i];  val[2] = 1/TA[i];
716:     col[3]     = net_start + 2*gbus[i]; col[4] = net_start + 2*gbus[i]+1;
717:     val[3]     = KA[i]*dVm_dVr/TA[i];         val[4] = KA[i]*dVm_dVi/TA[i];
718:     MatSetValues(J,1,row,5,col,val,INSERT_VALUES);
719:     idx        = idx + 9;
720:   }


723:   for (i=0; i<nbus; i++) {
724:     MatGetRow(user->Ybus,2*i,&ncols,&cols,&yvals);
725:     row[0] = net_start + 2*i;
726:     for (k=0; k<ncols; k++) {
727:       col[k] = net_start + cols[k];
728:       val[k] = yvals[k];
729:     }
730:     MatSetValues(J,1,row,ncols,col,val,INSERT_VALUES);
731:     MatRestoreRow(user->Ybus,2*i,&ncols,&cols,&yvals);

733:     MatGetRow(user->Ybus,2*i+1,&ncols,&cols,&yvals);
734:     row[0] = net_start + 2*i+1;
735:     for (k=0; k<ncols; k++) {
736:       col[k] = net_start + cols[k];
737:       val[k] = yvals[k];
738:     }
739:     MatSetValues(J,1,row,ncols,col,val,INSERT_VALUES);
740:     MatRestoreRow(user->Ybus,2*i+1,&ncols,&cols,&yvals);
741:   }

743:   MatAssemblyBegin(J,MAT_FLUSH_ASSEMBLY);
744:   MatAssemblyEnd(J,MAT_FLUSH_ASSEMBLY);


747:   VecGetArray(user->V0,&v0);
748:   for (i=0; i < nload; i++) {
749:     Vr      = xnet[2*lbus[i]]; /* Real part of load bus voltage */
750:     Vi      = xnet[2*lbus[i]+1]; /* Imaginary part of the load bus voltage */
751:     Vm      = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2= Vm*Vm; Vm4 = Vm2*Vm2;
752:     Vm0     = PetscSqrtScalar(v0[2*lbus[i]]*v0[2*lbus[i]] + v0[2*lbus[i]+1]*v0[2*lbus[i]+1]);
753:     PD      = QD = 0.0;
754:     dPD_dVr = dPD_dVi = dQD_dVr = dQD_dVi = 0.0;
755:     for (k=0; k < ld_nsegsp[i]; k++) {
756:       PD      += ld_alphap[k]*PD0[i]*PetscPowScalar((Vm/Vm0),ld_betap[k]);
757:       dPD_dVr += ld_alphap[k]*ld_betap[k]*PD0[i]*PetscPowScalar((1/Vm0),ld_betap[k])*Vr*PetscPowScalar(Vm,(ld_betap[k]-2));
758:       dPD_dVi += ld_alphap[k]*ld_betap[k]*PD0[i]*PetscPowScalar((1/Vm0),ld_betap[k])*Vi*PetscPowScalar(Vm,(ld_betap[k]-2));
759:     }
760:     for (k=0; k < ld_nsegsq[i]; k++) {
761:       QD      += ld_alphaq[k]*QD0[i]*PetscPowScalar((Vm/Vm0),ld_betaq[k]);
762:       dQD_dVr += ld_alphaq[k]*ld_betaq[k]*QD0[i]*PetscPowScalar((1/Vm0),ld_betaq[k])*Vr*PetscPowScalar(Vm,(ld_betaq[k]-2));
763:       dQD_dVi += ld_alphaq[k]*ld_betaq[k]*QD0[i]*PetscPowScalar((1/Vm0),ld_betaq[k])*Vi*PetscPowScalar(Vm,(ld_betaq[k]-2));
764:     }

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

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

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


776:     /*    fnet[2*lbus[i]]   += IDi; */
777:     row[0] = net_start + 2*lbus[i];
778:     col[0] = net_start + 2*lbus[i];  col[1] = net_start + 2*lbus[i]+1;
779:     val[0] = dIDi_dVr;               val[1] = dIDi_dVi;
780:     MatSetValues(J,1,row,2,col,val,ADD_VALUES);
781:     /*    fnet[2*lbus[i]+1] += IDr; */
782:     row[0] = net_start + 2*lbus[i]+1;
783:     col[0] = net_start + 2*lbus[i];  col[1] = net_start + 2*lbus[i]+1;
784:     val[0] = dIDr_dVr;               val[1] = dIDr_dVi;
785:     MatSetValues(J,1,row,2,col,val,ADD_VALUES);
786:   }
787:   VecRestoreArray(user->V0,&v0);

789:   VecRestoreArray(Xgen,&xgen);
790:   VecRestoreArray(Xnet,&xnet);

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

794:   MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
795:   MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
796:   return(0);
797: }

799: /*
800:    J = [I, 0
801:         dg_dx, dg_dy]
802: */
803: PetscErrorCode AlgJacobian(SNES snes,Vec X,Mat A,Mat B,void *ctx)
804: {
806:   Userctx        *user=(Userctx*)ctx;

809:   ResidualJacobian(snes,X,A,B,ctx);
810:   MatSetOption(A,MAT_KEEP_NONZERO_PATTERN,PETSC_TRUE);
811:   MatZeroRowsIS(A,user->is_diff,1.0,NULL,NULL);
812:   return(0);
813: }

815: /*
816:    J = [a*I-df_dx, -df_dy
817:         dg_dx, dg_dy]
818: */

820: PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,Userctx *user)
821: {
823:   SNES           snes;
824:   PetscScalar    atmp = (PetscScalar) a;
825:   PetscInt       i,row;

828:   user->t = t;

830:   TSGetSNES(ts,&snes);
831:   ResidualJacobian(snes,X,A,B,user);
832:   for (i=0;i < ngen;i++) {
833:     row = 9*i;
834:     MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
835:     row  = 9*i+1;
836:     MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
837:     row  = 9*i+2;
838:     MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
839:     row  = 9*i+3;
840:     MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
841:     row  = 9*i+6;
842:     MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
843:     row  = 9*i+7;
844:     MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
845:     row  = 9*i+8;
846:     MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
847:   }
848:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
849:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
850:   return(0);
851: }

853: /* Matrix JacobianP is constant so that it only needs to be evaluated once */
854: static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A, void *ctx0)
855: {
857:   PetscScalar    a;
858:   PetscInt       row,col;
859:   Userctx        *ctx=(Userctx*)ctx0;


863:   if (ctx->jacp_flg) {
864:     MatZeroEntries(A);

866:     for (col=0;col<3;col++) {
867:       a    = 1.0/M[col];
868:       row  = 9*col+3;
869:       MatSetValues(A,1,&row,1,&col,&a,INSERT_VALUES);
870:     }

872:     ctx->jacp_flg = PETSC_FALSE;

874:     MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
875:     MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
876:   }
877:   return(0);
878: }

880: static PetscErrorCode CostIntegrand(TS ts,PetscReal t,Vec U,Vec R,Userctx *user)
881: {
883:   PetscScalar    *u,*r;
884:   PetscInt       idx;
885:   Vec            Xgen,Xnet;
886:   PetscScalar    *xgen;
887:   PetscInt       i;

890:   DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
891:   DMCompositeScatter(user->dmpgrid,U,Xgen,Xnet);

893:   VecGetArray(Xgen,&xgen);

895:   VecGetArray(U,&u);
896:   VecGetArray(R,&r);
897:   r[0] = 0.;
898:   idx = 0;
899:   for (i=0;i<ngen;i++) {
900:     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);
901:     idx  += 9;
902:   }
903:   VecRestoreArray(R,&r);
904:   DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);
905:   return(0);
906: }

908: static PetscErrorCode DRDYFunction(TS ts,PetscReal t,Vec U,Vec *drdy,Userctx *user)
909: {
911:   Vec            Xgen,Xnet,Dgen,Dnet;
912:   PetscScalar    *xgen,*dgen;
913:   PetscInt       i;
914:   PetscInt       idx;

917:   DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
918:   DMCompositeGetLocalVectors(user->dmpgrid,&Dgen,&Dnet);
919:   DMCompositeScatter(user->dmpgrid,U,Xgen,Xnet);
920:   DMCompositeScatter(user->dmpgrid,drdy[0],Dgen,Dnet);

922:   VecGetArray(Xgen,&xgen);
923:   VecGetArray(Dgen,&dgen);

925:   idx = 0;
926:   for (i=0;i<ngen;i++) {
927:     dgen[idx+3] = 0.;
928:     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);
929:     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);
930:     idx += 9;
931:   }

933:   VecRestoreArray(Dgen,&dgen);
934:   DMCompositeGather(user->dmpgrid,INSERT_VALUES,drdy[0],Dgen,Dnet);
935:   DMCompositeRestoreLocalVectors(user->dmpgrid,&Dgen,&Dnet);
936:   DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);
937:   return(0);
938: }

940: static PetscErrorCode DRDPFunction(TS ts,PetscReal t,Vec U,Vec *drdp,Userctx *user)
941: {
943:   return(0);
944: }

946: PetscErrorCode ComputeSensiP(Vec lambda,Vec mu,Vec *DICDP,Userctx *user)
947: {
949:   PetscScalar    *x,*y,sensip;
950:   PetscInt       i;

953:   VecGetArray(lambda,&x);
954:   VecGetArray(mu,&y);

956:   for (i=0;i<3;i++) {
957:     VecDot(lambda,DICDP[i],&sensip);
958:     sensip = sensip+y[i];
959:     /* PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt %D th parameter: %g \n",i,(double)sensip); */
960:      y[i] = sensip;
961:   }
962:   VecRestoreArray(mu,&y);
963:   return(0);
964: }

966: int main(int argc,char **argv)
967: {
968:   Userctx            user;
969:   Vec                p;
970:   PetscScalar        *x_ptr;
971:   PetscErrorCode     ierr;
972:   PetscMPIInt        size;
973:   PetscInt           i;
974:   PetscViewer        Xview,Ybusview;
975:   PetscInt           *idx2;
976:   Tao                tao;
977:   KSP                ksp;
978:   PC                 pc;
979:   Vec                lowerb,upperb;

981:   PetscInitialize(&argc,&argv,"petscoptions",help);
982:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
983:   if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs");

985:   user.jacp_flg   = PETSC_TRUE;
986:   user.neqs_gen   = 9*ngen; /* # eqs. for generator subsystem */
987:   user.neqs_net   = 2*nbus; /* # eqs. for network subsystem   */
988:   user.neqs_pgrid = user.neqs_gen + user.neqs_net;

990:   /* Create indices for differential and algebraic equations */
991:   PetscMalloc1(7*ngen,&idx2);
992:   for (i=0; i<ngen; i++) {
993:     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;
994:     idx2[7*i+4] = 9*i+6; idx2[7*i+5] = 9*i+7; idx2[7*i+6] = 9*i+8;
995:   }
996:   ISCreateGeneral(PETSC_COMM_WORLD,7*ngen,idx2,PETSC_COPY_VALUES,&user.is_diff);
997:   ISComplement(user.is_diff,0,user.neqs_pgrid,&user.is_alg);
998:   PetscFree(idx2);

1000:   /* Set run time options */
1001:   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Transient stability fault options","");
1002:   {
1003:     user.tfaulton  = 1.0;
1004:     user.tfaultoff = 1.2;
1005:     user.Rfault    = 0.0001;
1006:     user.faultbus  = 8;
1007:     PetscOptionsReal("-tfaulton","","",user.tfaulton,&user.tfaulton,NULL);
1008:     PetscOptionsReal("-tfaultoff","","",user.tfaultoff,&user.tfaultoff,NULL);
1009:     PetscOptionsInt("-faultbus","","",user.faultbus,&user.faultbus,NULL);
1010:     user.t0        = 0.0;
1011:     user.tmax      = 1.3;
1012:     PetscOptionsReal("-t0","","",user.t0,&user.t0,NULL);
1013:     PetscOptionsReal("-tmax","","",user.tmax,&user.tmax,NULL);
1014:     user.freq_u    = 61.0;
1015:     user.freq_l    = 59.0;
1016:     user.pow       = 2;
1017:     PetscOptionsReal("-frequ","","",user.freq_u,&user.freq_u,NULL);
1018:     PetscOptionsReal("-freql","","",user.freq_l,&user.freq_l,NULL);
1019:     PetscOptionsInt("-pow","","",user.pow,&user.pow,NULL);

1021:   }
1022:   PetscOptionsEnd();

1024:   /* Create DMs for generator and network subsystems */
1025:   DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,user.neqs_gen,1,1,NULL,&user.dmgen);
1026:   DMSetOptionsPrefix(user.dmgen,"dmgen_");
1027:   DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,user.neqs_net,1,1,NULL,&user.dmnet);
1028:   DMSetOptionsPrefix(user.dmnet,"dmnet_");
1029:   DMSetFromOptions(user.dmnet);
1030:   DMSetUp(user.dmnet);
1031:   /* Create a composite DM packer and add the two DMs */
1032:   DMCompositeCreate(PETSC_COMM_WORLD,&user.dmpgrid);
1033:   DMSetOptionsPrefix(user.dmpgrid,"pgrid_");
1034:   DMSetFromOptions(user.dmgen);
1035:   DMSetUp(user.dmgen);
1036:   DMCompositeAddDM(user.dmpgrid,user.dmgen);
1037:   DMCompositeAddDM(user.dmpgrid,user.dmnet);

1039:   /* Read initial voltage vector and Ybus */
1040:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"X.bin",FILE_MODE_READ,&Xview);
1041:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"Ybus.bin",FILE_MODE_READ,&Ybusview);

1043:   VecCreate(PETSC_COMM_WORLD,&user.V0);
1044:   VecSetSizes(user.V0,PETSC_DECIDE,user.neqs_net);
1045:   VecLoad(user.V0,Xview);

1047:   MatCreate(PETSC_COMM_WORLD,&user.Ybus);
1048:   MatSetSizes(user.Ybus,PETSC_DECIDE,PETSC_DECIDE,user.neqs_net,user.neqs_net);
1049:   MatSetType(user.Ybus,MATBAIJ);
1050:   /*  MatSetBlockSize(ctx->Ybus,2); */
1051:   MatLoad(user.Ybus,Ybusview);

1053:   PetscViewerDestroy(&Xview);
1054:   PetscViewerDestroy(&Ybusview);

1056:   /* Allocate space for Jacobians */
1057:   MatCreate(PETSC_COMM_WORLD,&user.J);
1058:   MatSetSizes(user.J,PETSC_DECIDE,PETSC_DECIDE,user.neqs_pgrid,user.neqs_pgrid);
1059:   MatSetFromOptions(user.J);
1060:   PreallocateJacobian(user.J,&user);

1062:   MatCreate(PETSC_COMM_WORLD,&user.Jacp);
1063:   MatSetSizes(user.Jacp,PETSC_DECIDE,PETSC_DECIDE,user.neqs_pgrid,3);
1064:   MatSetFromOptions(user.Jacp);
1065:   MatSetUp(user.Jacp);
1066:   MatZeroEntries(user.Jacp); /* initialize to zeros */

1068:   /* Create TAO solver and set desired solution method */
1069:   TaoCreate(PETSC_COMM_WORLD,&tao);
1070:   TaoSetType(tao,TAOBLMVM);
1071:   /*
1072:      Optimization starts
1073:   */
1074:   /* Set initial solution guess */
1075:   VecCreateSeq(PETSC_COMM_WORLD,3,&p);
1076:   VecGetArray(p,&x_ptr);
1077:   x_ptr[0] = PG[0]; x_ptr[1] = PG[1]; x_ptr[2] = PG[2];
1078:   VecRestoreArray(p,&x_ptr);

1080:   TaoSetInitialVector(tao,p);
1081:   /* Set routine for function and gradient evaluation */
1082:   TaoSetObjectiveAndGradientRoutine(tao,FormFunctionGradient,&user);

1084:   /* Set bounds for the optimization */
1085:   VecDuplicate(p,&lowerb);
1086:   VecDuplicate(p,&upperb);
1087:   VecGetArray(lowerb,&x_ptr);
1088:   x_ptr[0] = 0.5; x_ptr[1] = 0.5; x_ptr[2] = 0.5;
1089:   VecRestoreArray(lowerb,&x_ptr);
1090:   VecGetArray(upperb,&x_ptr);
1091:   x_ptr[0] = 2.0; x_ptr[1] = 2.0; x_ptr[2] = 2.0;
1092:   VecRestoreArray(upperb,&x_ptr);
1093:   TaoSetVariableBounds(tao,lowerb,upperb);

1095:   /* Check for any TAO command line options */
1096:   TaoSetFromOptions(tao);
1097:   TaoGetKSP(tao,&ksp);
1098:   if (ksp) {
1099:     KSPGetPC(ksp,&pc);
1100:     PCSetType(pc,PCNONE);
1101:   }

1103:   /* SOLVE THE APPLICATION */
1104:   TaoSolve(tao);

1106:   VecView(p,PETSC_VIEWER_STDOUT_WORLD);
1107:   /* Free TAO data structures */
1108:   TaoDestroy(&tao);

1110:   DMDestroy(&user.dmgen);
1111:   DMDestroy(&user.dmnet);
1112:   DMDestroy(&user.dmpgrid);
1113:   ISDestroy(&user.is_diff);
1114:   ISDestroy(&user.is_alg);

1116:   MatDestroy(&user.J);
1117:   MatDestroy(&user.Jacp);
1118:   MatDestroy(&user.Ybus);
1119:   /* MatDestroy(&user.Sol); */
1120:   VecDestroy(&user.V0);
1121:   VecDestroy(&p);
1122:   VecDestroy(&lowerb);
1123:   VecDestroy(&upperb);
1124:   PetscFinalize();
1125:   return ierr;
1126: }

1128: /* ------------------------------------------------------------------ */
1129: /*
1130:    FormFunction - Evaluates the function and corresponding gradient.

1132:    Input Parameters:
1133:    tao - the Tao context
1134:    X   - the input vector
1135:    ptr - optional user-defined context, as set by TaoSetObjectiveAndGradientRoutine()

1137:    Output Parameters:
1138:    f   - the newly evaluated function
1139:    G   - the newly evaluated gradient
1140: */
1141: PetscErrorCode FormFunctionGradient(Tao tao,Vec P,PetscReal *f,Vec G,void *ctx0)
1142: {
1143:   TS             ts;
1144:   SNES           snes_alg;
1146:   Userctx        *ctx = (Userctx*)ctx0;
1147:   Vec            X;
1148:   PetscInt       i;
1149:   /* sensitivity context */
1150:   PetscScalar    *x_ptr;
1151:   Vec            lambda[1],q;
1152:   Vec            mu[1];
1153:   PetscInt       steps1,steps2,steps3;
1154:   Vec            DICDP[3];
1155:   Vec            F_alg;
1156:   PetscInt       row_loc,col_loc;
1157:   PetscScalar    val;
1158:   Vec            Xdot;

1161:   VecGetArrayRead(P,(const PetscScalar**)&x_ptr);
1162:   PG[0] = x_ptr[0];
1163:   PG[1] = x_ptr[1];
1164:   PG[2] = x_ptr[2];
1165:   VecRestoreArrayRead(P,(const PetscScalar**)&x_ptr);

1167:   ctx->stepnum = 0;

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

1171:   /* Create matrix to save solutions at each time step */
1172:   /* MatCreateSeqDense(PETSC_COMM_SELF,ctx->neqs_pgrid+1,1002,NULL,&ctx->Sol); */
1173:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1174:      Create timestepping solver context
1175:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1176:   TSCreate(PETSC_COMM_WORLD,&ts);
1177:   TSSetProblemType(ts,TS_NONLINEAR);
1178:   TSSetType(ts,TSCN);
1179:   TSSetIFunction(ts,NULL,(TSIFunction) IFunction,ctx);
1180:   TSSetIJacobian(ts,ctx->J,ctx->J,(TSIJacobian)IJacobian,ctx);
1181:   TSSetApplicationContext(ts,ctx);

1183:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1184:      Set initial conditions
1185:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1186:   SetInitialGuess(X,ctx);

1188:   /* Approximate DICDP with finite difference, we want to zero out network variables */
1189:   for (i=0;i<3;i++) {
1190:     VecDuplicate(X,&DICDP[i]);
1191:   }
1192:   DICDPFiniteDifference(X,DICDP,ctx);

1194:   VecDuplicate(X,&F_alg);
1195:   SNESCreate(PETSC_COMM_WORLD,&snes_alg);
1196:   SNESSetFunction(snes_alg,F_alg,AlgFunction,ctx);
1197:   MatZeroEntries(ctx->J);
1198:   SNESSetJacobian(snes_alg,ctx->J,ctx->J,AlgJacobian,ctx);
1199:   SNESSetOptionsPrefix(snes_alg,"alg_");
1200:   SNESSetFromOptions(snes_alg);
1201:   ctx->alg_flg = PETSC_TRUE;
1202:   /* Solve the algebraic equations */
1203:   SNESSolve(snes_alg,NULL,X);

1205:   /* Just to set up the Jacobian structure */
1206:   VecDuplicate(X,&Xdot);
1207:   IJacobian(ts,0.0,X,Xdot,0.0,ctx->J,ctx->J,ctx);
1208:   VecDestroy(&Xdot);

1210:   ctx->stepnum++;

1212:   /*
1213:     Save trajectory of solution so that TSAdjointSolve() may be used
1214:   */
1215:   TSSetSaveTrajectory(ts);

1217:   TSSetTimeStep(ts,0.01);
1218:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);
1219:   TSSetFromOptions(ts);
1220:   /* TSSetPostStep(ts,SaveSolution); */


1223:   /* Prefault period */
1224:   ctx->alg_flg = PETSC_FALSE;
1225:   TSSetTime(ts,0.0);
1226:   TSSetMaxTime(ts,ctx->tfaulton);
1227:   TSSolve(ts,X);
1228:   TSGetStepNumber(ts,&steps1);

1230:   /* Create the nonlinear solver for solving the algebraic system */
1231:   /* Note that although the algebraic system needs to be solved only for
1232:      Idq and V, we reuse the entire system including xgen. The xgen
1233:      variables are held constant by setting their residuals to 0 and
1234:      putting a 1 on the Jacobian diagonal for xgen rows
1235:   */
1236:   MatZeroEntries(ctx->J);

1238:   /* Apply disturbance - resistive fault at ctx->faultbus */
1239:   /* This is done by adding shunt conductance to the diagonal location
1240:      in the Ybus matrix */
1241:   row_loc = 2*ctx->faultbus; col_loc = 2*ctx->faultbus+1; /* Location for G */
1242:   val     = 1/ctx->Rfault;
1243:   MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);
1244:   row_loc = 2*ctx->faultbus+1; col_loc = 2*ctx->faultbus; /* Location for G */
1245:   val     = 1/ctx->Rfault;
1246:   MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);

1248:   MatAssemblyBegin(ctx->Ybus,MAT_FINAL_ASSEMBLY);
1249:   MatAssemblyEnd(ctx->Ybus,MAT_FINAL_ASSEMBLY);

1251:   ctx->alg_flg = PETSC_TRUE;
1252:   /* Solve the algebraic equations */
1253:   SNESSolve(snes_alg,NULL,X);

1255:   ctx->stepnum++;

1257:   /* Disturbance period */
1258:   ctx->alg_flg = PETSC_FALSE;
1259:   TSSetTime(ts,ctx->tfaulton);
1260:   TSSetMaxTime(ts,ctx->tfaultoff);
1261:   TSSolve(ts,X);
1262:   TSGetStepNumber(ts,&steps2);
1263:   steps2 -= steps1;

1265:   /* Remove the fault */
1266:   row_loc = 2*ctx->faultbus; col_loc = 2*ctx->faultbus+1;
1267:   val     = -1/ctx->Rfault;
1268:   MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);
1269:   row_loc = 2*ctx->faultbus+1; col_loc = 2*ctx->faultbus;
1270:   val     = -1/ctx->Rfault;
1271:   MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);

1273:   MatAssemblyBegin(ctx->Ybus,MAT_FINAL_ASSEMBLY);
1274:   MatAssemblyEnd(ctx->Ybus,MAT_FINAL_ASSEMBLY);

1276:   MatZeroEntries(ctx->J);

1278:   ctx->alg_flg = PETSC_TRUE;

1280:   /* Solve the algebraic equations */
1281:   SNESSolve(snes_alg,NULL,X);

1283:   ctx->stepnum++;

1285:   /* Post-disturbance period */
1286:   ctx->alg_flg = PETSC_TRUE;
1287:   TSSetTime(ts,ctx->tfaultoff);
1288:   TSSetMaxTime(ts,ctx->tmax);
1289:   TSSolve(ts,X);
1290:   TSGetStepNumber(ts,&steps3);
1291:   steps3 -= steps2;
1292:   steps3 -= steps1;

1294:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1295:      Adjoint model starts here
1296:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1297:   TSSetPostStep(ts,NULL);
1298:   MatCreateVecs(ctx->J,&lambda[0],NULL);
1299:   /*   Set initial conditions for the adjoint integration */
1300:   VecZeroEntries(lambda[0]);

1302:   MatCreateVecs(ctx->Jacp,&mu[0],NULL);
1303:   VecZeroEntries(mu[0]);
1304:   TSSetCostGradients(ts,1,lambda,mu);

1306:   /*   Set RHS JacobianP */
1307:   TSSetRHSJacobianP(ts,ctx->Jacp,RHSJacobianP,ctx);

1309:   TSSetCostIntegrand(ts,1,NULL,(PetscErrorCode (*)(TS,PetscReal,Vec,Vec,void*))CostIntegrand,
1310:                                         (PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*))DRDYFunction,
1311:                                         (PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*))DRDPFunction,PETSC_FALSE,ctx);

1313:   TSAdjointSetSteps(ts,steps3);
1314:   TSAdjointSolve(ts);

1316:   MatZeroEntries(ctx->J);
1317:   /* Applying disturbance - resistive fault at ctx->faultbus */
1318:   /* This is done by deducting shunt conductance to the diagonal location
1319:      in the Ybus matrix */
1320:   row_loc = 2*ctx->faultbus; col_loc = 2*ctx->faultbus+1; /* Location for G */
1321:   val     = 1./ctx->Rfault;
1322:   MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);
1323:   row_loc = 2*ctx->faultbus+1; col_loc = 2*ctx->faultbus; /* Location for G */
1324:   val     = 1./ctx->Rfault;
1325:   MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);

1327:   MatAssemblyBegin(ctx->Ybus,MAT_FINAL_ASSEMBLY);
1328:   MatAssemblyEnd(ctx->Ybus,MAT_FINAL_ASSEMBLY);


1331:   /*   Set number of steps for the adjoint integration */
1332:   TSAdjointSetSteps(ts,steps2);
1333:   TSAdjointSolve(ts);

1335:   MatZeroEntries(ctx->J);
1336:   /* remove the fault */
1337:   row_loc = 2*ctx->faultbus; col_loc = 2*ctx->faultbus+1; /* Location for G */
1338:   val     = -1./ctx->Rfault;
1339:   MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);
1340:   row_loc = 2*ctx->faultbus+1; col_loc = 2*ctx->faultbus; /* Location for G */
1341:   val     = -1./ctx->Rfault;
1342:   MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);

1344:   MatAssemblyBegin(ctx->Ybus,MAT_FINAL_ASSEMBLY);
1345:   MatAssemblyEnd(ctx->Ybus,MAT_FINAL_ASSEMBLY);

1347:   /*   Set number of steps for the adjoint integration */
1348:   TSAdjointSetSteps(ts,steps1);
1349:   TSAdjointSolve(ts);


1352:   ComputeSensiP(lambda[0],mu[0],DICDP,ctx);
1353:   VecCopy(mu[0],G);
1354:   TSGetCostIntegral(ts,&q);
1355:   VecGetArray(q,&x_ptr);
1356:   *f   = x_ptr[0];

1358:   VecDestroy(&lambda[0]);
1359:   VecDestroy(&mu[0]);

1361:   SNESDestroy(&snes_alg);
1362:   VecDestroy(&F_alg);
1363:   VecDestroy(&X);
1364:   TSDestroy(&ts);
1365:   for (i=0;i<3;i++) {
1366:     VecDestroy(&DICDP[i]);
1367:   }
1368:   return(0);
1369: }

1371: /*TEST

1373:    build:
1374:       requires: double !complex !define(PETSC_USE_64BIT_INDICES)

1376:    test:
1377:       args: -viewer_binary_skip_info -tao_monitor -tao_gttol .2
1378:       localrunfiles: petscoptions X.bin Ybus.bin

1380: TEST*/