Actual source code: ex9busopt.c

petsc-3.12.0 2019-09-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:   Mat         DRDU,DRDP;
108: } Userctx;


111: /* Converts from machine frame (dq) to network (phase a real,imag) reference frame */
112: PetscErrorCode dq2ri(PetscScalar Fd,PetscScalar Fq,PetscScalar delta,PetscScalar *Fr, PetscScalar *Fi)
113: {
115:   *Fr =  Fd*PetscSinScalar(delta) + Fq*PetscCosScalar(delta);
116:   *Fi = -Fd*PetscCosScalar(delta) + Fq*PetscSinScalar(delta);
117:   return(0);
118: }

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

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

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

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

168:   M[0] = 2*H[0]/w_s; M[1] = 2*H[1]/w_s; M[2] = 2*H[2]/w_s;
169:   D[0] = 0.1*M[0]; D[1] = 0.1*M[1]; D[2] = 0.1*M[2];

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

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

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

180:   for (i=0; i < ngen; i++) {
181:     Vr  = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */
182:     Vi  = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */
183:     Vm  = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm;
184:     IGr = (Vr*PG[i] + Vi*QG[i])/Vm2;
185:     IGi = (Vi*PG[i] - Vr*QG[i])/Vm2;

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

189:     theta = PETSC_PI/2.0 - delta;

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

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

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

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

202:     /* The generator variables are ordered as [Eqp,Edp,delta,w,Id,Iq] */
203:     xgen[idx]   = Eqp;
204:     xgen[idx+1] = Edp;
205:     xgen[idx+2] = delta;
206:     xgen[idx+3] = w_s;

208:     idx = idx + 4;

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

213:     idx = idx + 2;

215:     /* Exciter */
216:     Efd = Eqp + (Xd[i] - Xdp[i])*Id;
217:     SE  = k1[i]*PetscExpScalar(k2[i]*Efd);
218:     VR  =  KE[i]*Efd + SE;
219:     RF  =  KF[i]*Efd/TF[i];

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

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

227:     idx = idx + 3;
228:   }

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

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

239: PetscErrorCode InitialGuess(Vec X,Userctx *user, const PetscScalar PGv[])
240: {
242:   Vec            Xgen,Xnet;
243:   PetscScalar    *xgen,*xnet;
244:   PetscInt       i,idx=0;
245:   PetscScalar    Vr,Vi,IGr,IGi,Vm,Vm2;
246:   PetscScalar    Eqp,Edp,delta;
247:   PetscScalar    Efd,RF,VR; /* Exciter variables */
248:   PetscScalar    Id,Iq;  /* Generator dq axis currents */
249:   PetscScalar    theta,Vd,Vq,SE;

252:   M[0] = 2*H[0]/w_s; M[1] = 2*H[1]/w_s; M[2] = 2*H[2]/w_s;
253:   D[0] = 0.1*M[0]; D[1] = 0.1*M[1]; D[2] = 0.1*M[2];

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

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

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

264:   for (i=0; i < ngen; i++) {
265:     Vr  = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */
266:     Vi  = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */
267:     Vm  = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm;
268:     IGr = (Vr*PGv[i] + Vi*QG[i])/Vm2;
269:     IGi = (Vi*PGv[i] - Vr*QG[i])/Vm2;

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

273:     theta = PETSC_PI/2.0 - delta;

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

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

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

284:     /* The generator variables are ordered as [Eqp,Edp,delta,w,Id,Iq] */
285:     xgen[idx]   = Eqp;
286:     xgen[idx+1] = Edp;
287:     xgen[idx+2] = delta;
288:     xgen[idx+3] = w_s;

290:     idx = idx + 4;

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

295:     idx = idx + 2;

297:     /* Exciter */
298:     Efd = Eqp + (Xd[i] - Xdp[i])*Id;
299:     SE  = k1[i]*PetscExpScalar(k2[i]*Efd);
300:     VR  =  KE[i]*Efd + SE;
301:     RF  =  KF[i]*Efd/TF[i];

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

307:     idx = idx + 3;
308:   }

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

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

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

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

329:   for (i=0;i<ngen;i++) {
330:     for (j=0;j<3;j++) PGv[j] = PG[j];
331:     PGv[i] = PG[i]+eps;
332:     InitialGuess(Y,user,PGv);
333:     InitialGuess(X,user,PG);

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


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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

431:     idx = idx + 9;
432:   }

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

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

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

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

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

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

476:   user->t = t;

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

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

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

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

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

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

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

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

551:     idx = idx + 9;
552:   }

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

657:     /*    fgen[idx+4] = Zdq_inv[0]*(-Edp + Vd) + Zdq_inv[1]*(-Eqp + Vq) + Id; */
658:     row[0] = idx+4;
659:     col[0] = idx;         col[1] = idx+1;        col[2] = idx + 2;
660:     val[0] = -Zdq_inv[1]; val[1] = -Zdq_inv[0];  val[2] = Zdq_inv[0]*dVd_ddelta + Zdq_inv[1]*dVq_ddelta;
661:     col[3] = idx + 4; col[4] = net_start+2*gbus[i];                     col[5] = net_start + 2*gbus[i]+1;
662:     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;
663:     MatSetValues(J,1,row,6,col,val,INSERT_VALUES);

665:     /*  fgen[idx+5] = Zdq_inv[2]*(-Edp + Vd) + Zdq_inv[3]*(-Eqp + Vq) + Iq; */
666:     row[0] = idx+5;
667:     col[0] = idx;         col[1] = idx+1;        col[2] = idx + 2;
668:     val[0] = -Zdq_inv[3]; val[1] = -Zdq_inv[2];  val[2] = Zdq_inv[2]*dVd_ddelta + Zdq_inv[3]*dVq_ddelta;
669:     col[3] = idx + 5; col[4] = net_start+2*gbus[i];                     col[5] = net_start + 2*gbus[i]+1;
670:     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;
671:     MatSetValues(J,1,row,6,col,val,INSERT_VALUES);

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

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

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

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

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

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

701:     /* Exciter differential equations */

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

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


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

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

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


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

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

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

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


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

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

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

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

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

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

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

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

829:   user->t = t;

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

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


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

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

873:     ctx->jacp_flg = PETSC_FALSE;

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

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

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

894:   VecGetArray(Xgen,&xgen);

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

909: static PetscErrorCode DRDUJacobianTranspose(TS ts,PetscReal t,Vec U,Mat DRDU,Mat B,Userctx *user)
910: {
912:   Vec            Xgen,Xnet,Dgen,Dnet;
913:   PetscScalar    *xgen,*dgen;
914:   PetscInt       i;
915:   PetscInt       idx;
916:   Vec            drdu_col;
917:   PetscScalar    *xarr;

920:   VecDuplicate(U,&drdu_col);
921:   MatDenseGetColumn(DRDU,0,&xarr);
922:   VecPlaceArray(drdu_col,xarr);
923:   DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
924:   DMCompositeGetLocalVectors(user->dmpgrid,&Dgen,&Dnet);
925:   DMCompositeScatter(user->dmpgrid,U,Xgen,Xnet);
926:   DMCompositeScatter(user->dmpgrid,drdu_col,Dgen,Dnet);

928:   VecGetArray(Xgen,&xgen);
929:   VecGetArray(Dgen,&dgen);

931:   idx = 0;
932:   for (i=0;i<ngen;i++) {
933:     dgen[idx+3] = 0.;
934:     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);
935:     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);
936:     idx += 9;
937:   }

939:   VecRestoreArray(Dgen,&dgen);
940:   VecRestoreArray(Xgen,&xgen);
941:   DMCompositeGather(user->dmpgrid,INSERT_VALUES,drdu_col,Dgen,Dnet);
942:   DMCompositeRestoreLocalVectors(user->dmpgrid,&Dgen,&Dnet);
943:   DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);
944:   VecResetArray(drdu_col);
945:   MatDenseRestoreColumn(DRDU,&xarr);
946:   VecDestroy(&drdu_col);
947:   return(0);
948: }

950: static PetscErrorCode DRDPJacobianTranspose(TS ts,PetscReal t,Vec U,Mat drdp,Userctx *user)
951: {
953:   return(0);
954: }

956: PetscErrorCode ComputeSensiP(Vec lambda,Vec mu,Vec *DICDP,Userctx *user)
957: {
959:   PetscScalar    *x,*y,sensip;
960:   PetscInt       i;

963:   VecGetArray(lambda,&x);
964:   VecGetArray(mu,&y);

966:   for (i=0;i<3;i++) {
967:     VecDot(lambda,DICDP[i],&sensip);
968:     sensip = sensip+y[i];
969:     /* PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt %D th parameter: %g \n",i,(double)sensip); */
970:      y[i] = sensip;
971:   }
972:   VecRestoreArray(mu,&y);
973:   return(0);
974: }

976: int main(int argc,char **argv)
977: {
978:   Userctx            user;
979:   Vec                p;
980:   PetscScalar        *x_ptr;
981:   PetscErrorCode     ierr;
982:   PetscMPIInt        size;
983:   PetscInt           i;
984:   PetscViewer        Xview,Ybusview;
985:   PetscInt           *idx2;
986:   Tao                tao;
987:   KSP                ksp;
988:   PC                 pc;
989:   Vec                lowerb,upperb;

991:   PetscInitialize(&argc,&argv,"petscoptions",help);if (ierr) return ierr;
992:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
993:   if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs");

995:   user.jacp_flg   = PETSC_TRUE;
996:   user.neqs_gen   = 9*ngen; /* # eqs. for generator subsystem */
997:   user.neqs_net   = 2*nbus; /* # eqs. for network subsystem   */
998:   user.neqs_pgrid = user.neqs_gen + user.neqs_net;

1000:   /* Create indices for differential and algebraic equations */
1001:   PetscMalloc1(7*ngen,&idx2);
1002:   for (i=0; i<ngen; i++) {
1003:     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;
1004:     idx2[7*i+4] = 9*i+6; idx2[7*i+5] = 9*i+7; idx2[7*i+6] = 9*i+8;
1005:   }
1006:   ISCreateGeneral(PETSC_COMM_WORLD,7*ngen,idx2,PETSC_COPY_VALUES,&user.is_diff);
1007:   ISComplement(user.is_diff,0,user.neqs_pgrid,&user.is_alg);
1008:   PetscFree(idx2);

1010:   /* Set run time options */
1011:   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Transient stability fault options","");
1012:   {
1013:     user.tfaulton  = 1.0;
1014:     user.tfaultoff = 1.2;
1015:     user.Rfault    = 0.0001;
1016:     user.faultbus  = 8;
1017:     PetscOptionsReal("-tfaulton","","",user.tfaulton,&user.tfaulton,NULL);
1018:     PetscOptionsReal("-tfaultoff","","",user.tfaultoff,&user.tfaultoff,NULL);
1019:     PetscOptionsInt("-faultbus","","",user.faultbus,&user.faultbus,NULL);
1020:     user.t0        = 0.0;
1021:     user.tmax      = 1.3;
1022:     PetscOptionsReal("-t0","","",user.t0,&user.t0,NULL);
1023:     PetscOptionsReal("-tmax","","",user.tmax,&user.tmax,NULL);
1024:     user.freq_u    = 61.0;
1025:     user.freq_l    = 59.0;
1026:     user.pow       = 2;
1027:     PetscOptionsReal("-frequ","","",user.freq_u,&user.freq_u,NULL);
1028:     PetscOptionsReal("-freql","","",user.freq_l,&user.freq_l,NULL);
1029:     PetscOptionsInt("-pow","","",user.pow,&user.pow,NULL);

1031:   }
1032:   PetscOptionsEnd();

1034:   /* Create DMs for generator and network subsystems */
1035:   DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,user.neqs_gen,1,1,NULL,&user.dmgen);
1036:   DMSetOptionsPrefix(user.dmgen,"dmgen_");
1037:   DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,user.neqs_net,1,1,NULL,&user.dmnet);
1038:   DMSetOptionsPrefix(user.dmnet,"dmnet_");
1039:   DMSetFromOptions(user.dmnet);
1040:   DMSetUp(user.dmnet);
1041:   /* Create a composite DM packer and add the two DMs */
1042:   DMCompositeCreate(PETSC_COMM_WORLD,&user.dmpgrid);
1043:   DMSetOptionsPrefix(user.dmpgrid,"pgrid_");
1044:   DMSetFromOptions(user.dmgen);
1045:   DMSetUp(user.dmgen);
1046:   DMCompositeAddDM(user.dmpgrid,user.dmgen);
1047:   DMCompositeAddDM(user.dmpgrid,user.dmnet);

1049:   /* Read initial voltage vector and Ybus */
1050:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"X.bin",FILE_MODE_READ,&Xview);
1051:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"Ybus.bin",FILE_MODE_READ,&Ybusview);

1053:   VecCreate(PETSC_COMM_WORLD,&user.V0);
1054:   VecSetSizes(user.V0,PETSC_DECIDE,user.neqs_net);
1055:   VecLoad(user.V0,Xview);

1057:   MatCreate(PETSC_COMM_WORLD,&user.Ybus);
1058:   MatSetSizes(user.Ybus,PETSC_DECIDE,PETSC_DECIDE,user.neqs_net,user.neqs_net);
1059:   MatSetType(user.Ybus,MATBAIJ);
1060:   /*  MatSetBlockSize(ctx->Ybus,2); */
1061:   MatLoad(user.Ybus,Ybusview);

1063:   PetscViewerDestroy(&Xview);
1064:   PetscViewerDestroy(&Ybusview);

1066:   /* Allocate space for Jacobians */
1067:   MatCreate(PETSC_COMM_WORLD,&user.J);
1068:   MatSetSizes(user.J,PETSC_DECIDE,PETSC_DECIDE,user.neqs_pgrid,user.neqs_pgrid);
1069:   MatSetFromOptions(user.J);
1070:   PreallocateJacobian(user.J,&user);

1072:   MatCreate(PETSC_COMM_WORLD,&user.Jacp);
1073:   MatSetSizes(user.Jacp,PETSC_DECIDE,PETSC_DECIDE,user.neqs_pgrid,3);
1074:   MatSetFromOptions(user.Jacp);
1075:   MatSetUp(user.Jacp);
1076:   MatZeroEntries(user.Jacp); /* initialize to zeros */

1078:   MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,3,1,NULL,&user.DRDP);
1079:   MatSetUp(user.DRDP);
1080:   MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,user.neqs_pgrid,1,NULL,&user.DRDU);
1081:   MatSetUp(user.DRDU);

1083:   /* Create TAO solver and set desired solution method */
1084:   TaoCreate(PETSC_COMM_WORLD,&tao);
1085:   TaoSetType(tao,TAOBLMVM);
1086:   /*
1087:      Optimization starts
1088:   */
1089:   /* Set initial solution guess */
1090:   VecCreateSeq(PETSC_COMM_WORLD,3,&p);
1091:   VecGetArray(p,&x_ptr);
1092:   x_ptr[0] = PG[0]; x_ptr[1] = PG[1]; x_ptr[2] = PG[2];
1093:   VecRestoreArray(p,&x_ptr);

1095:   TaoSetInitialVector(tao,p);
1096:   /* Set routine for function and gradient evaluation */
1097:   TaoSetObjectiveAndGradientRoutine(tao,FormFunctionGradient,&user);

1099:   /* Set bounds for the optimization */
1100:   VecDuplicate(p,&lowerb);
1101:   VecDuplicate(p,&upperb);
1102:   VecGetArray(lowerb,&x_ptr);
1103:   x_ptr[0] = 0.5; x_ptr[1] = 0.5; x_ptr[2] = 0.5;
1104:   VecRestoreArray(lowerb,&x_ptr);
1105:   VecGetArray(upperb,&x_ptr);
1106:   x_ptr[0] = 2.0; x_ptr[1] = 2.0; x_ptr[2] = 2.0;
1107:   VecRestoreArray(upperb,&x_ptr);
1108:   TaoSetVariableBounds(tao,lowerb,upperb);

1110:   /* Check for any TAO command line options */
1111:   TaoSetFromOptions(tao);
1112:   TaoGetKSP(tao,&ksp);
1113:   if (ksp) {
1114:     KSPGetPC(ksp,&pc);
1115:     PCSetType(pc,PCNONE);
1116:   }

1118:   /* SOLVE THE APPLICATION */
1119:   TaoSolve(tao);

1121:   VecView(p,PETSC_VIEWER_STDOUT_WORLD);
1122:   /* Free TAO data structures */
1123:   TaoDestroy(&tao);

1125:   DMDestroy(&user.dmgen);
1126:   DMDestroy(&user.dmnet);
1127:   DMDestroy(&user.dmpgrid);
1128:   ISDestroy(&user.is_diff);
1129:   ISDestroy(&user.is_alg);

1131:   MatDestroy(&user.J);
1132:   MatDestroy(&user.Jacp);
1133:   MatDestroy(&user.Ybus);
1134:   /* MatDestroy(&user.Sol); */
1135:   VecDestroy(&user.V0);
1136:   VecDestroy(&p);
1137:   VecDestroy(&lowerb);
1138:   VecDestroy(&upperb);
1139:   MatDestroy(&user.DRDU);
1140:   MatDestroy(&user.DRDP);
1141:   PetscFinalize();
1142:   return ierr;
1143: }

1145: /* ------------------------------------------------------------------ */
1146: /*
1147:    FormFunction - Evaluates the function and corresponding gradient.

1149:    Input Parameters:
1150:    tao - the Tao context
1151:    X   - the input vector
1152:    ptr - optional user-defined context, as set by TaoSetObjectiveAndGradientRoutine()

1154:    Output Parameters:
1155:    f   - the newly evaluated function
1156:    G   - the newly evaluated gradient
1157: */
1158: PetscErrorCode FormFunctionGradient(Tao tao,Vec P,PetscReal *f,Vec G,void *ctx0)
1159: {
1160:   TS             ts,quadts;
1161:   SNES           snes_alg;
1163:   Userctx        *ctx = (Userctx*)ctx0;
1164:   Vec            X;
1165:   PetscInt       i;
1166:   /* sensitivity context */
1167:   PetscScalar    *x_ptr;
1168:   Vec            lambda[1],q;
1169:   Vec            mu[1];
1170:   PetscInt       steps1,steps2,steps3;
1171:   Vec            DICDP[3];
1172:   Vec            F_alg;
1173:   PetscInt       row_loc,col_loc;
1174:   PetscScalar    val;
1175:   Vec            Xdot;

1178:   VecGetArrayRead(P,(const PetscScalar**)&x_ptr);
1179:   PG[0] = x_ptr[0];
1180:   PG[1] = x_ptr[1];
1181:   PG[2] = x_ptr[2];
1182:   VecRestoreArrayRead(P,(const PetscScalar**)&x_ptr);

1184:   ctx->stepnum = 0;

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

1188:   /* Create matrix to save solutions at each time step */
1189:   /* MatCreateSeqDense(PETSC_COMM_SELF,ctx->neqs_pgrid+1,1002,NULL,&ctx->Sol); */
1190:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1191:      Create timestepping solver context
1192:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1193:   TSCreate(PETSC_COMM_WORLD,&ts);
1194:   TSSetProblemType(ts,TS_NONLINEAR);
1195:   TSSetType(ts,TSCN);
1196:   TSSetIFunction(ts,NULL,(TSIFunction) IFunction,ctx);
1197:   TSSetIJacobian(ts,ctx->J,ctx->J,(TSIJacobian)IJacobian,ctx);
1198:   TSSetApplicationContext(ts,ctx);
1199:   /*   Set RHS JacobianP */
1200:   TSSetRHSJacobianP(ts,ctx->Jacp,RHSJacobianP,ctx);

1202:   TSCreateQuadratureTS(ts,PETSC_FALSE,&quadts);
1203:   TSSetRHSFunction(quadts,NULL,(TSRHSFunction)CostIntegrand,ctx);
1204:   TSSetRHSJacobian(quadts,ctx->DRDU,ctx->DRDU,(TSRHSJacobian)DRDUJacobianTranspose,ctx);
1205:   TSSetRHSJacobianP(quadts,ctx->DRDP,(TSRHSJacobianP)DRDPJacobianTranspose,ctx);

1207:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1208:      Set initial conditions
1209:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1210:   SetInitialGuess(X,ctx);

1212:   /* Approximate DICDP with finite difference, we want to zero out network variables */
1213:   for (i=0;i<3;i++) {
1214:     VecDuplicate(X,&DICDP[i]);
1215:   }
1216:   DICDPFiniteDifference(X,DICDP,ctx);

1218:   VecDuplicate(X,&F_alg);
1219:   SNESCreate(PETSC_COMM_WORLD,&snes_alg);
1220:   SNESSetFunction(snes_alg,F_alg,AlgFunction,ctx);
1221:   MatZeroEntries(ctx->J);
1222:   SNESSetJacobian(snes_alg,ctx->J,ctx->J,AlgJacobian,ctx);
1223:   SNESSetOptionsPrefix(snes_alg,"alg_");
1224:   SNESSetFromOptions(snes_alg);
1225:   ctx->alg_flg = PETSC_TRUE;
1226:   /* Solve the algebraic equations */
1227:   SNESSolve(snes_alg,NULL,X);

1229:   /* Just to set up the Jacobian structure */
1230:   VecDuplicate(X,&Xdot);
1231:   IJacobian(ts,0.0,X,Xdot,0.0,ctx->J,ctx->J,ctx);
1232:   VecDestroy(&Xdot);

1234:   ctx->stepnum++;

1236:   /*
1237:     Save trajectory of solution so that TSAdjointSolve() may be used
1238:   */
1239:   TSSetSaveTrajectory(ts);

1241:   TSSetTimeStep(ts,0.01);
1242:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);
1243:   TSSetFromOptions(ts);
1244:   /* TSSetPostStep(ts,SaveSolution); */


1247:   /* Prefault period */
1248:   ctx->alg_flg = PETSC_FALSE;
1249:   TSSetTime(ts,0.0);
1250:   TSSetMaxTime(ts,ctx->tfaulton);
1251:   TSSolve(ts,X);
1252:   TSGetStepNumber(ts,&steps1);

1254:   /* Create the nonlinear solver for solving the algebraic system */
1255:   /* Note that although the algebraic system needs to be solved only for
1256:      Idq and V, we reuse the entire system including xgen. The xgen
1257:      variables are held constant by setting their residuals to 0 and
1258:      putting a 1 on the Jacobian diagonal for xgen rows
1259:   */
1260:   MatZeroEntries(ctx->J);

1262:   /* Apply disturbance - resistive fault at ctx->faultbus */
1263:   /* This is done by adding shunt conductance to the diagonal location
1264:      in the Ybus matrix */
1265:   row_loc = 2*ctx->faultbus; col_loc = 2*ctx->faultbus+1; /* Location for G */
1266:   val     = 1/ctx->Rfault;
1267:   MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);
1268:   row_loc = 2*ctx->faultbus+1; col_loc = 2*ctx->faultbus; /* Location for G */
1269:   val     = 1/ctx->Rfault;
1270:   MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);

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

1275:   ctx->alg_flg = PETSC_TRUE;
1276:   /* Solve the algebraic equations */
1277:   SNESSolve(snes_alg,NULL,X);

1279:   ctx->stepnum++;

1281:   /* Disturbance period */
1282:   ctx->alg_flg = PETSC_FALSE;
1283:   TSSetTime(ts,ctx->tfaulton);
1284:   TSSetMaxTime(ts,ctx->tfaultoff);
1285:   TSSolve(ts,X);
1286:   TSGetStepNumber(ts,&steps2);
1287:   steps2 -= steps1;

1289:   /* Remove the fault */
1290:   row_loc = 2*ctx->faultbus; col_loc = 2*ctx->faultbus+1;
1291:   val     = -1/ctx->Rfault;
1292:   MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);
1293:   row_loc = 2*ctx->faultbus+1; col_loc = 2*ctx->faultbus;
1294:   val     = -1/ctx->Rfault;
1295:   MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);

1297:   MatAssemblyBegin(ctx->Ybus,MAT_FINAL_ASSEMBLY);
1298:   MatAssemblyEnd(ctx->Ybus,MAT_FINAL_ASSEMBLY);

1300:   MatZeroEntries(ctx->J);

1302:   ctx->alg_flg = PETSC_TRUE;

1304:   /* Solve the algebraic equations */
1305:   SNESSolve(snes_alg,NULL,X);

1307:   ctx->stepnum++;

1309:   /* Post-disturbance period */
1310:   ctx->alg_flg = PETSC_TRUE;
1311:   TSSetTime(ts,ctx->tfaultoff);
1312:   TSSetMaxTime(ts,ctx->tmax);
1313:   TSSolve(ts,X);
1314:   TSGetStepNumber(ts,&steps3);
1315:   steps3 -= steps2;
1316:   steps3 -= steps1;

1318:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1319:      Adjoint model starts here
1320:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1321:   TSSetPostStep(ts,NULL);
1322:   MatCreateVecs(ctx->J,&lambda[0],NULL);
1323:   /*   Set initial conditions for the adjoint integration */
1324:   VecZeroEntries(lambda[0]);

1326:   MatCreateVecs(ctx->Jacp,&mu[0],NULL);
1327:   VecZeroEntries(mu[0]);
1328:   TSSetCostGradients(ts,1,lambda,mu);

1330:   TSAdjointSetSteps(ts,steps3);
1331:   TSAdjointSolve(ts);

1333:   MatZeroEntries(ctx->J);
1334:   /* Applying disturbance - resistive fault at ctx->faultbus */
1335:   /* This is done by deducting shunt conductance to the diagonal location
1336:      in the Ybus matrix */
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);


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

1352:   MatZeroEntries(ctx->J);
1353:   /* remove the fault */
1354:   row_loc = 2*ctx->faultbus; col_loc = 2*ctx->faultbus+1; /* Location for G */
1355:   val     = -1./ctx->Rfault;
1356:   MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);
1357:   row_loc = 2*ctx->faultbus+1; col_loc = 2*ctx->faultbus; /* Location for G */
1358:   val     = -1./ctx->Rfault;
1359:   MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);

1361:   MatAssemblyBegin(ctx->Ybus,MAT_FINAL_ASSEMBLY);
1362:   MatAssemblyEnd(ctx->Ybus,MAT_FINAL_ASSEMBLY);

1364:   /*   Set number of steps for the adjoint integration */
1365:   TSAdjointSetSteps(ts,steps1);
1366:   TSAdjointSolve(ts);

1368:   ComputeSensiP(lambda[0],mu[0],DICDP,ctx);
1369:   VecCopy(mu[0],G);

1371:   TSGetQuadratureTS(ts,NULL,&quadts);
1372:   TSGetSolution(quadts,&q);
1373:   VecGetArray(q,&x_ptr);
1374:   *f   = x_ptr[0];
1375:   x_ptr[0] = 0;
1376:   VecRestoreArray(q,&x_ptr);

1378:   VecDestroy(&lambda[0]);
1379:   VecDestroy(&mu[0]);

1381:   SNESDestroy(&snes_alg);
1382:   VecDestroy(&F_alg);
1383:   VecDestroy(&X);
1384:   TSDestroy(&ts);
1385:   for (i=0;i<3;i++) {
1386:     VecDestroy(&DICDP[i]);
1387:   }
1388:   return(0);
1389: }

1391: /*TEST

1393:    build:
1394:       requires: double !complex !define(PETSC_USE_64BIT_INDICES)

1396:    test:
1397:       args: -viewer_binary_skip_info -tao_monitor -tao_gttol .2
1398:       localrunfiles: petscoptions X.bin Ybus.bin

1400: TEST*/