Actual source code: ex9busoptfd.c

  1: static char help[] = "Using finite difference for the problem in ex9busopt.c \n\n";

  3: /*
  4:   Use finite difference approximations to solve the same optimization problem as in ex9busopt.c.
  5:  */

  7: #include <petsctao.h>
  8: #include <petscts.h>
  9: #include <petscdm.h>
 10: #include <petscdmda.h>
 11: #include <petscdmcomposite.h>

 13: PetscErrorCode FormFunction(Tao,Vec,PetscReal*,void*);

 15: #define freq 60
 16: #define w_s (2*PETSC_PI*freq)

 18: /* Sizes and indices */
 19: const PetscInt nbus    = 9; /* Number of network buses */
 20: const PetscInt ngen    = 3; /* Number of generators */
 21: const PetscInt nload   = 3; /* Number of loads */
 22: const PetscInt gbus[3] = {0,1,2}; /* Buses at which generators are incident */
 23: const PetscInt lbus[3] = {4,5,7}; /* Buses at which loads are incident */

 25: /* Generator real and reactive powers (found via loadflow) */
 26: PetscScalar PG[3] = { 0.69,1.59,0.69};
 27: /* PetscScalar PG[3] = {0.716786142395021,1.630000000000000,0.850000000000000};*/
 28: const PetscScalar QG[3] = {0.270702180178785,0.066120127797275,-0.108402221791588};
 29: /* Generator constants */
 30: const PetscScalar H[3]    = {23.64,6.4,3.01};   /* Inertia constant */
 31: const PetscScalar Rs[3]   = {0.0,0.0,0.0}; /* Stator Resistance */
 32: const PetscScalar Xd[3]   = {0.146,0.8958,1.3125};  /* d-axis reactance */
 33: const PetscScalar Xdp[3]  = {0.0608,0.1198,0.1813}; /* d-axis transient reactance */
 34: 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 */
 35: const PetscScalar Xqp[3]  = {0.0969,0.1969,0.25}; /* q-axis transient reactance */
 36: const PetscScalar Td0p[3] = {8.96,6.0,5.89}; /* d-axis open circuit time constant */
 37: const PetscScalar Tq0p[3] = {0.31,0.535,0.6}; /* q-axis open circuit time constant */
 38: PetscScalar M[3]; /* M = 2*H/w_s */
 39: PetscScalar D[3]; /* D = 0.1*M */

 41: PetscScalar TM[3]; /* Mechanical Torque */
 42: /* Exciter system constants */
 43: const PetscScalar KA[3] = {20.0,20.0,20.0};  /* Voltage regulartor gain constant */
 44: const PetscScalar TA[3] = {0.2,0.2,0.2};     /* Voltage regulator time constant */
 45: const PetscScalar KE[3] = {1.0,1.0,1.0};     /* Exciter gain constant */
 46: const PetscScalar TE[3] = {0.314,0.314,0.314}; /* Exciter time constant */
 47: const PetscScalar KF[3] = {0.063,0.063,0.063};  /* Feedback stabilizer gain constant */
 48: const PetscScalar TF[3] = {0.35,0.35,0.35};    /* Feedback stabilizer time constant */
 49: const PetscScalar k1[3] = {0.0039,0.0039,0.0039};
 50: const PetscScalar k2[3] = {1.555,1.555,1.555};  /* k1 and k2 for calculating the saturation function SE = k1*exp(k2*Efd) */

 52: PetscScalar Vref[3];
 53: /* Load constants
 54:   We use a composite load model that describes the load and reactive powers at each time instant as follows
 55:   P(t) = \sum\limits_{i=0}^ld_nsegsp \ld_alphap_i*P_D0(\frac{V_m(t)}{V_m0})^\ld_betap_i
 56:   Q(t) = \sum\limits_{i=0}^ld_nsegsq \ld_alphaq_i*Q_D0(\frac{V_m(t)}{V_m0})^\ld_betaq_i
 57:   where
 58:     ld_nsegsp,ld_nsegsq - Number of individual load models for real and reactive power loads
 59:     ld_alphap,ld_alphap - Percentage contribution (weights) or loads
 60:     P_D0                - Real power load
 61:     Q_D0                - Reactive power load
 62:     V_m(t)              - Voltage magnitude at time t
 63:     V_m0                - Voltage magnitude at t = 0
 64:     ld_betap, ld_betaq  - exponents describing the load model for real and reactive part

 66:     Note: All loads have the same characteristic currently.
 67: */
 68: const PetscScalar PD0[3] = {1.25,0.9,1.0};
 69: const PetscScalar QD0[3] = {0.5,0.3,0.35};
 70: const PetscInt    ld_nsegsp[3] = {3,3,3};
 71: const PetscScalar ld_alphap[3] = {1.0,0.0,0.0};
 72: const PetscScalar ld_betap[3]  = {2.0,1.0,0.0};
 73: const PetscInt    ld_nsegsq[3] = {3,3,3};
 74: const PetscScalar ld_alphaq[3] = {1.0,0.0,0.0};
 75: const PetscScalar ld_betaq[3]  = {2.0,1.0,0.0};

 77: typedef struct {
 78:   DM          dmgen, dmnet; /* DMs to manage generator and network subsystem */
 79:   DM          dmpgrid; /* Composite DM to manage the entire power grid */
 80:   Mat         Ybus; /* Network admittance matrix */
 81:   Vec         V0;  /* Initial voltage vector (Power flow solution) */
 82:   PetscReal   tfaulton,tfaultoff; /* Fault on and off times */
 83:   PetscInt    faultbus; /* Fault bus */
 84:   PetscScalar Rfault;
 85:   PetscReal   t0,tmax;
 86:   PetscInt    neqs_gen,neqs_net,neqs_pgrid;
 87:   Mat         Sol; /* Matrix to save solution at each time step */
 88:   PetscInt    stepnum;
 89:   PetscBool   alg_flg;
 90:   PetscReal   t;
 91:   IS          is_diff; /* indices for differential equations */
 92:   IS          is_alg; /* indices for algebraic equations */
 93:   PetscReal   freq_u,freq_l; /* upper and lower frequency limit */
 94:   PetscInt    pow; /* power coefficient used in the cost function */
 95:   Vec         vec_q;
 96: } Userctx;

 98: /* Converts from machine frame (dq) to network (phase a real,imag) reference frame */
 99: PetscErrorCode dq2ri(PetscScalar Fd,PetscScalar Fq,PetscScalar delta,PetscScalar *Fr, PetscScalar *Fi)
100: {
102:   *Fr =  Fd*PetscSinScalar(delta) + Fq*PetscCosScalar(delta);
103:   *Fi = -Fd*PetscCosScalar(delta) + Fq*PetscSinScalar(delta);
104:   return(0);
105: }

107: /* Converts from network frame ([phase a real,imag) to machine (dq) reference frame */
108: PetscErrorCode ri2dq(PetscScalar Fr,PetscScalar Fi,PetscScalar delta,PetscScalar *Fd, PetscScalar *Fq)
109: {
111:   *Fd =  Fr*PetscSinScalar(delta) - Fi*PetscCosScalar(delta);
112:   *Fq =  Fr*PetscCosScalar(delta) + Fi*PetscSinScalar(delta);
113:   return(0);
114: }

116: PetscErrorCode SetInitialGuess(Vec X,Userctx *user)
117: {
119:   Vec            Xgen,Xnet;
120:   PetscScalar    *xgen,*xnet;
121:   PetscInt       i,idx=0;
122:   PetscScalar    Vr,Vi,IGr,IGi,Vm,Vm2;
123:   PetscScalar    Eqp,Edp,delta;
124:   PetscScalar    Efd,RF,VR; /* Exciter variables */
125:   PetscScalar    Id,Iq;  /* Generator dq axis currents */
126:   PetscScalar    theta,Vd,Vq,SE;

129:   M[0] = 2*H[0]/w_s; M[1] = 2*H[1]/w_s; M[2] = 2*H[2]/w_s;
130:   D[0] = 0.1*M[0]; D[1] = 0.1*M[1]; D[2] = 0.1*M[2];

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

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

137:   /* Generator subsystem initialization */
138:   VecGetArray(Xgen,&xgen);
139:   VecGetArray(Xnet,&xnet);

141:   for (i=0; i < ngen; i++) {
142:     Vr  = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */
143:     Vi  = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */
144:     Vm  = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm;
145:     IGr = (Vr*PG[i] + Vi*QG[i])/Vm2;
146:     IGi = (Vi*PG[i] - Vr*QG[i])/Vm2;

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

150:     theta = PETSC_PI/2.0 - delta;

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

155:     Vd = Vr*PetscCosScalar(theta) - Vi*PetscSinScalar(theta);
156:     Vq = Vr*PetscSinScalar(theta) + Vi*PetscCosScalar(theta);

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

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

163:     /* The generator variables are ordered as [Eqp,Edp,delta,w,Id,Iq] */
164:     xgen[idx]   = Eqp;
165:     xgen[idx+1] = Edp;
166:     xgen[idx+2] = delta;
167:     xgen[idx+3] = w_s;

169:     idx = idx + 4;

171:     xgen[idx]   = Id;
172:     xgen[idx+1] = Iq;

174:     idx = idx + 2;

176:     /* Exciter */
177:     Efd = Eqp + (Xd[i] - Xdp[i])*Id;
178:     SE  = k1[i]*PetscExpScalar(k2[i]*Efd);
179:     VR  =  KE[i]*Efd + SE;
180:     RF  =  KF[i]*Efd/TF[i];

182:     xgen[idx]   = Efd;
183:     xgen[idx+1] = RF;
184:     xgen[idx+2] = VR;

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

188:     idx = idx + 3;
189:   }

191:   VecRestoreArray(Xgen,&xgen);
192:   VecRestoreArray(Xnet,&xnet);

194:   /* VecView(Xgen,0); */
195:   DMCompositeGather(user->dmpgrid,INSERT_VALUES,X,Xgen,Xnet);
196:   DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);
197:   return(0);
198: }

200: /* Computes F = [-f(x,y);g(x,y)] */
201: PetscErrorCode ResidualFunction(SNES snes,Vec X, Vec F, Userctx *user)
202: {
204:   Vec            Xgen,Xnet,Fgen,Fnet;
205:   PetscScalar    *xgen,*xnet,*fgen,*fnet;
206:   PetscInt       i,idx=0;
207:   PetscScalar    Vr,Vi,Vm,Vm2;
208:   PetscScalar    Eqp,Edp,delta,w; /* Generator variables */
209:   PetscScalar    Efd,RF,VR; /* Exciter variables */
210:   PetscScalar    Id,Iq;  /* Generator dq axis currents */
211:   PetscScalar    Vd,Vq,SE;
212:   PetscScalar    IGr,IGi,IDr,IDi;
213:   PetscScalar    Zdq_inv[4],det;
214:   PetscScalar    PD,QD,Vm0,*v0;
215:   PetscInt       k;

218:   VecZeroEntries(F);
219:   DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
220:   DMCompositeGetLocalVectors(user->dmpgrid,&Fgen,&Fnet);
221:   DMCompositeScatter(user->dmpgrid,X,Xgen,Xnet);
222:   DMCompositeScatter(user->dmpgrid,F,Fgen,Fnet);

224:   /* Network current balance residual IG + Y*V + IL = 0. Only YV is added here.
225:      The generator current injection, IG, and load current injection, ID are added later
226:   */
227:   /* Note that the values in Ybus are stored assuming the imaginary current balance
228:      equation is ordered first followed by real current balance equation for each bus.
229:      Thus imaginary current contribution goes in location 2*i, and
230:      real current contribution in 2*i+1
231:   */
232:   MatMult(user->Ybus,Xnet,Fnet);

234:   VecGetArray(Xgen,&xgen);
235:   VecGetArray(Xnet,&xnet);
236:   VecGetArray(Fgen,&fgen);
237:   VecGetArray(Fnet,&fnet);

239:   /* Generator subsystem */
240:   for (i=0; i < ngen; i++) {
241:     Eqp   = xgen[idx];
242:     Edp   = xgen[idx+1];
243:     delta = xgen[idx+2];
244:     w     = xgen[idx+3];
245:     Id    = xgen[idx+4];
246:     Iq    = xgen[idx+5];
247:     Efd   = xgen[idx+6];
248:     RF    = xgen[idx+7];
249:     VR    = xgen[idx+8];

251:     /* Generator differential equations */
252:     fgen[idx]   = (Eqp + (Xd[i] - Xdp[i])*Id - Efd)/Td0p[i];
253:     fgen[idx+1] = (Edp - (Xq[i] - Xqp[i])*Iq)/Tq0p[i];
254:     fgen[idx+2] = -w + w_s;
255:     fgen[idx+3] = (-TM[i] + Edp*Id + Eqp*Iq + (Xqp[i] - Xdp[i])*Id*Iq + D[i]*(w - w_s))/M[i];

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

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

264:     Zdq_inv[0] = Rs[i]/det;
265:     Zdq_inv[1] = Xqp[i]/det;
266:     Zdq_inv[2] = -Xdp[i]/det;
267:     Zdq_inv[3] = Rs[i]/det;

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

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

275:     fnet[2*gbus[i]]   -= IGi;
276:     fnet[2*gbus[i]+1] -= IGr;

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

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

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

287:     idx = idx + 9;
288:   }

290:   VecGetArray(user->V0,&v0);
291:   for (i=0; i < nload; i++) {
292:     Vr  = xnet[2*lbus[i]]; /* Real part of load bus voltage */
293:     Vi  = xnet[2*lbus[i]+1]; /* Imaginary part of the load bus voltage */
294:     Vm  = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm;
295:     Vm0 = PetscSqrtScalar(v0[2*lbus[i]]*v0[2*lbus[i]] + v0[2*lbus[i]+1]*v0[2*lbus[i]+1]);
296:     PD  = QD = 0.0;
297:     for (k=0; k < ld_nsegsp[i]; k++) PD += ld_alphap[k]*PD0[i]*PetscPowScalar((Vm/Vm0),ld_betap[k]);
298:     for (k=0; k < ld_nsegsq[i]; k++) QD += ld_alphaq[k]*QD0[i]*PetscPowScalar((Vm/Vm0),ld_betaq[k]);

300:     /* Load currents */
301:     IDr = (PD*Vr + QD*Vi)/Vm2;
302:     IDi = (-QD*Vr + PD*Vi)/Vm2;

304:     fnet[2*lbus[i]]   += IDi;
305:     fnet[2*lbus[i]+1] += IDr;
306:   }
307:   VecRestoreArray(user->V0,&v0);

309:   VecRestoreArray(Xgen,&xgen);
310:   VecRestoreArray(Xnet,&xnet);
311:   VecRestoreArray(Fgen,&fgen);
312:   VecRestoreArray(Fnet,&fnet);

314:   DMCompositeGather(user->dmpgrid,INSERT_VALUES,F,Fgen,Fnet);
315:   DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);
316:   DMCompositeRestoreLocalVectors(user->dmpgrid,&Fgen,&Fnet);
317:   return(0);
318: }

320: /* \dot{x} - f(x,y)
321:      g(x,y) = 0
322:  */
323: PetscErrorCode IFunction(TS ts,PetscReal t, Vec X, Vec Xdot, Vec F, Userctx *user)
324: {
325:   PetscErrorCode    ierr;
326:   SNES              snes;
327:   PetscScalar       *f;
328:   const PetscScalar *xdot;
329:   PetscInt          i;

332:   user->t = t;

334:   TSGetSNES(ts,&snes);
335:   ResidualFunction(snes,X,F,user);
336:   VecGetArray(F,&f);
337:   VecGetArrayRead(Xdot,&xdot);
338:   for (i=0;i < ngen;i++) {
339:     f[9*i]   += xdot[9*i];
340:     f[9*i+1] += xdot[9*i+1];
341:     f[9*i+2] += xdot[9*i+2];
342:     f[9*i+3] += xdot[9*i+3];
343:     f[9*i+6] += xdot[9*i+6];
344:     f[9*i+7] += xdot[9*i+7];
345:     f[9*i+8] += xdot[9*i+8];
346:   }
347:   VecRestoreArray(F,&f);
348:   VecRestoreArrayRead(Xdot,&xdot);
349:   return(0);
350: }

352: /* This function is used for solving the algebraic system only during fault on and
353:    off times. It computes the entire F and then zeros out the part corresponding to
354:    differential equations
355:  F = [0;g(y)];
356: */
357: PetscErrorCode AlgFunction(SNES snes, Vec X, Vec F, void *ctx)
358: {
360:   Userctx        *user=(Userctx*)ctx;
361:   PetscScalar    *f;
362:   PetscInt       i;

365:   ResidualFunction(snes,X,F,user);
366:   VecGetArray(F,&f);
367:   for (i=0; i < ngen; i++) {
368:     f[9*i]   = 0;
369:     f[9*i+1] = 0;
370:     f[9*i+2] = 0;
371:     f[9*i+3] = 0;
372:     f[9*i+6] = 0;
373:     f[9*i+7] = 0;
374:     f[9*i+8] = 0;
375:   }
376:   VecRestoreArray(F,&f);
377:   return(0);
378: }

380: PetscErrorCode PreallocateJacobian(Mat J, Userctx *user)
381: {
383:   PetscInt       *d_nnz;
384:   PetscInt       i,idx=0,start=0;
385:   PetscInt       ncols;

388:   PetscMalloc1(user->neqs_pgrid,&d_nnz);
389:   for (i=0; i<user->neqs_pgrid; i++) d_nnz[i] = 0;
390:   /* Generator subsystem */
391:   for (i=0; i < ngen; i++) {

393:     d_nnz[idx]   += 3;
394:     d_nnz[idx+1] += 2;
395:     d_nnz[idx+2] += 2;
396:     d_nnz[idx+3] += 5;
397:     d_nnz[idx+4] += 6;
398:     d_nnz[idx+5] += 6;

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

403:     d_nnz[idx+6] += 2;
404:     d_nnz[idx+7] += 2;
405:     d_nnz[idx+8] += 5;

407:     idx = idx + 9;
408:   }

410:   start = user->neqs_gen;

412:   for (i=0; i < nbus; i++) {
413:     MatGetRow(user->Ybus,2*i,&ncols,NULL,NULL);
414:     d_nnz[start+2*i]   += ncols;
415:     d_nnz[start+2*i+1] += ncols;
416:     MatRestoreRow(user->Ybus,2*i,&ncols,NULL,NULL);
417:   }

419:   MatSeqAIJSetPreallocation(J,0,d_nnz);

421:   PetscFree(d_nnz);
422:   return(0);
423: }

425: /*
426:    J = [-df_dx, -df_dy
427:         dg_dx, dg_dy]
428: */
429: PetscErrorCode ResidualJacobian(SNES snes,Vec X,Mat J,Mat B,void *ctx)
430: {
431:   PetscErrorCode    ierr;
432:   Userctx           *user=(Userctx*)ctx;
433:   Vec               Xgen,Xnet;
434:   PetscScalar       *xgen,*xnet;
435:   PetscInt          i,idx=0;
436:   PetscScalar       Vr,Vi,Vm,Vm2;
437:   PetscScalar       Eqp,Edp,delta; /* Generator variables */
438:   PetscScalar       Efd; /* Exciter variables */
439:   PetscScalar       Id,Iq;  /* Generator dq axis currents */
440:   PetscScalar       Vd,Vq;
441:   PetscScalar       val[10];
442:   PetscInt          row[2],col[10];
443:   PetscInt          net_start=user->neqs_gen;
444:   PetscScalar       Zdq_inv[4],det;
445:   PetscScalar       dVd_dVr,dVd_dVi,dVq_dVr,dVq_dVi,dVd_ddelta,dVq_ddelta;
446:   PetscScalar       dIGr_ddelta,dIGi_ddelta,dIGr_dId,dIGr_dIq,dIGi_dId,dIGi_dIq;
447:   PetscScalar       dSE_dEfd;
448:   PetscScalar       dVm_dVd,dVm_dVq,dVm_dVr,dVm_dVi;
449:   PetscInt          ncols;
450:   const PetscInt    *cols;
451:   const PetscScalar *yvals;
452:   PetscInt          k;
453:   PetscScalar       PD,QD,Vm0,*v0,Vm4;
454:   PetscScalar       dPD_dVr,dPD_dVi,dQD_dVr,dQD_dVi;
455:   PetscScalar       dIDr_dVr,dIDr_dVi,dIDi_dVr,dIDi_dVi;

458:   MatZeroEntries(B);
459:   DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
460:   DMCompositeScatter(user->dmpgrid,X,Xgen,Xnet);

462:   VecGetArray(Xgen,&xgen);
463:   VecGetArray(Xnet,&xnet);

465:   /* Generator subsystem */
466:   for (i=0; i < ngen; i++) {
467:     Eqp   = xgen[idx];
468:     Edp   = xgen[idx+1];
469:     delta = xgen[idx+2];
470:     Id    = xgen[idx+4];
471:     Iq    = xgen[idx+5];
472:     Efd   = xgen[idx+6];

474:     /*    fgen[idx]   = (Eqp + (Xd[i] - Xdp[i])*Id - Efd)/Td0p[i]; */
475:     row[0] = idx;
476:     col[0] = idx;           col[1] = idx+4;          col[2] = idx+6;
477:     val[0] = 1/ Td0p[i]; val[1] = (Xd[i] - Xdp[i])/ Td0p[i]; val[2] = -1/Td0p[i];

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

481:     /*    fgen[idx+1] = (Edp - (Xq[i] - Xqp[i])*Iq)/Tq0p[i]; */
482:     row[0] = idx + 1;
483:     col[0] = idx + 1;       col[1] = idx+5;
484:     val[0] = 1/Tq0p[i]; val[1] = -(Xq[i] - Xqp[i])/Tq0p[i];
485:     MatSetValues(J,1,row,2,col,val,INSERT_VALUES);

487:     /*    fgen[idx+2] = - w + w_s; */
488:     row[0] = idx + 2;
489:     col[0] = idx + 2; col[1] = idx + 3;
490:     val[0] = 0;       val[1] = -1;
491:     MatSetValues(J,1,row,2,col,val,INSERT_VALUES);

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

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

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

505:     Zdq_inv[0] = Rs[i]/det;
506:     Zdq_inv[1] = Xqp[i]/det;
507:     Zdq_inv[2] = -Xdp[i]/det;
508:     Zdq_inv[3] = Rs[i]/det;

510:     dVd_dVr    = PetscSinScalar(delta); dVd_dVi = -PetscCosScalar(delta);
511:     dVq_dVr    = PetscCosScalar(delta); dVq_dVi = PetscSinScalar(delta);
512:     dVd_ddelta = Vr*PetscCosScalar(delta) + Vi*PetscSinScalar(delta);
513:     dVq_ddelta = -Vr*PetscSinScalar(delta) + Vi*PetscCosScalar(delta);

515:     /*    fgen[idx+4] = Zdq_inv[0]*(-Edp + Vd) + Zdq_inv[1]*(-Eqp + Vq) + Id; */
516:     row[0] = idx+4;
517:     col[0] = idx;         col[1] = idx+1;        col[2] = idx + 2;
518:     val[0] = -Zdq_inv[1]; val[1] = -Zdq_inv[0];  val[2] = Zdq_inv[0]*dVd_ddelta + Zdq_inv[1]*dVq_ddelta;
519:     col[3] = idx + 4; col[4] = net_start+2*gbus[i];                     col[5] = net_start + 2*gbus[i]+1;
520:     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;
521:     MatSetValues(J,1,row,6,col,val,INSERT_VALUES);

523:     /*  fgen[idx+5] = Zdq_inv[2]*(-Edp + Vd) + Zdq_inv[3]*(-Eqp + Vq) + Iq; */
524:     row[0] = idx+5;
525:     col[0] = idx;         col[1] = idx+1;        col[2] = idx + 2;
526:     val[0] = -Zdq_inv[3]; val[1] = -Zdq_inv[2];  val[2] = Zdq_inv[2]*dVd_ddelta + Zdq_inv[3]*dVq_ddelta;
527:     col[3] = idx + 5; col[4] = net_start+2*gbus[i];                     col[5] = net_start + 2*gbus[i]+1;
528:     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;
529:     MatSetValues(J,1,row,6,col,val,INSERT_VALUES);

531:     dIGr_ddelta = Id*PetscCosScalar(delta) - Iq*PetscSinScalar(delta);
532:     dIGi_ddelta = Id*PetscSinScalar(delta) + Iq*PetscCosScalar(delta);
533:     dIGr_dId    = PetscSinScalar(delta);  dIGr_dIq = PetscCosScalar(delta);
534:     dIGi_dId    = -PetscCosScalar(delta); dIGi_dIq = PetscSinScalar(delta);

536:     /* fnet[2*gbus[i]]   -= IGi; */
537:     row[0] = net_start + 2*gbus[i];
538:     col[0] = idx+2;        col[1] = idx + 4;   col[2] = idx + 5;
539:     val[0] = -dIGi_ddelta; val[1] = -dIGi_dId; val[2] = -dIGi_dIq;
540:     MatSetValues(J,1,row,3,col,val,INSERT_VALUES);

542:     /* fnet[2*gbus[i]+1]   -= IGr; */
543:     row[0] = net_start + 2*gbus[i]+1;
544:     col[0] = idx+2;        col[1] = idx + 4;   col[2] = idx + 5;
545:     val[0] = -dIGr_ddelta; val[1] = -dIGr_dId; val[2] = -dIGr_dIq;
546:     MatSetValues(J,1,row,3,col,val,INSERT_VALUES);

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

550:     /*    fgen[idx+6] = (KE[i]*Efd + SE - VR)/TE[i]; */
551:     /*    SE  = k1[i]*PetscExpScalar(k2[i]*Efd); */

553:     dSE_dEfd = k1[i]*k2[i]*PetscExpScalar(k2[i]*Efd);

555:     row[0] = idx + 6;
556:     col[0] = idx + 6;                     col[1] = idx + 8;
557:     val[0] = (KE[i] + dSE_dEfd)/TE[i];  val[1] = -1/TE[i];
558:     MatSetValues(J,1,row,2,col,val,INSERT_VALUES);

560:     /* Exciter differential equations */

562:     /*    fgen[idx+7] = (RF - KF[i]*Efd/TF[i])/TF[i]; */
563:     row[0] = idx + 7;
564:     col[0] = idx + 6;       col[1] = idx + 7;
565:     val[0] = (-KF[i]/TF[i])/TF[i];  val[1] = 1/TF[i];
566:     MatSetValues(J,1,row,2,col,val,INSERT_VALUES);

568:     /*    fgen[idx+8] = (VR - KA[i]*RF + KA[i]*KF[i]*Efd/TF[i] - KA[i]*(Vref[i] - Vm))/TA[i]; */
569:     /* Vm = (Vd^2 + Vq^2)^0.5; */
570:     dVm_dVd    = Vd/Vm; dVm_dVq = Vq/Vm;
571:     dVm_dVr    = dVm_dVd*dVd_dVr + dVm_dVq*dVq_dVr;
572:     dVm_dVi    = dVm_dVd*dVd_dVi + dVm_dVq*dVq_dVi;
573:     row[0]     = idx + 8;
574:     col[0]     = idx + 6;           col[1] = idx + 7; col[2] = idx + 8;
575:     val[0]     = (KA[i]*KF[i]/TF[i])/TA[i]; val[1] = -KA[i]/TA[i];  val[2] = 1/TA[i];
576:     col[3]     = net_start + 2*gbus[i]; col[4] = net_start + 2*gbus[i]+1;
577:     val[3]     = KA[i]*dVm_dVr/TA[i];         val[4] = KA[i]*dVm_dVi/TA[i];
578:     MatSetValues(J,1,row,5,col,val,INSERT_VALUES);
579:     idx        = idx + 9;
580:   }

582:   for (i=0; i<nbus; i++) {
583:     MatGetRow(user->Ybus,2*i,&ncols,&cols,&yvals);
584:     row[0] = net_start + 2*i;
585:     for (k=0; k<ncols; k++) {
586:       col[k] = net_start + cols[k];
587:       val[k] = yvals[k];
588:     }
589:     MatSetValues(J,1,row,ncols,col,val,INSERT_VALUES);
590:     MatRestoreRow(user->Ybus,2*i,&ncols,&cols,&yvals);

592:     MatGetRow(user->Ybus,2*i+1,&ncols,&cols,&yvals);
593:     row[0] = net_start + 2*i+1;
594:     for (k=0; k<ncols; k++) {
595:       col[k] = net_start + cols[k];
596:       val[k] = yvals[k];
597:     }
598:     MatSetValues(J,1,row,ncols,col,val,INSERT_VALUES);
599:     MatRestoreRow(user->Ybus,2*i+1,&ncols,&cols,&yvals);
600:   }

602:   MatAssemblyBegin(J,MAT_FLUSH_ASSEMBLY);
603:   MatAssemblyEnd(J,MAT_FLUSH_ASSEMBLY);

605:   VecGetArray(user->V0,&v0);
606:   for (i=0; i < nload; i++) {
607:     Vr      = xnet[2*lbus[i]]; /* Real part of load bus voltage */
608:     Vi      = xnet[2*lbus[i]+1]; /* Imaginary part of the load bus voltage */
609:     Vm      = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm; Vm4 = Vm2*Vm2;
610:     Vm0     = PetscSqrtScalar(v0[2*lbus[i]]*v0[2*lbus[i]] + v0[2*lbus[i]+1]*v0[2*lbus[i]+1]);
611:     PD      = QD = 0.0;
612:     dPD_dVr = dPD_dVi = dQD_dVr = dQD_dVi = 0.0;
613:     for (k=0; k < ld_nsegsp[i]; k++) {
614:       PD      += ld_alphap[k]*PD0[i]*PetscPowScalar((Vm/Vm0),ld_betap[k]);
615:       dPD_dVr += ld_alphap[k]*ld_betap[k]*PD0[i]*PetscPowScalar((1/Vm0),ld_betap[k])*Vr*PetscPowScalar(Vm,(ld_betap[k]-2));
616:       dPD_dVi += ld_alphap[k]*ld_betap[k]*PD0[i]*PetscPowScalar((1/Vm0),ld_betap[k])*Vi*PetscPowScalar(Vm,(ld_betap[k]-2));
617:     }
618:     for (k=0; k < ld_nsegsq[i]; k++) {
619:       QD      += ld_alphaq[k]*QD0[i]*PetscPowScalar((Vm/Vm0),ld_betaq[k]);
620:       dQD_dVr += ld_alphaq[k]*ld_betaq[k]*QD0[i]*PetscPowScalar((1/Vm0),ld_betaq[k])*Vr*PetscPowScalar(Vm,(ld_betaq[k]-2));
621:       dQD_dVi += ld_alphaq[k]*ld_betaq[k]*QD0[i]*PetscPowScalar((1/Vm0),ld_betaq[k])*Vi*PetscPowScalar(Vm,(ld_betaq[k]-2));
622:     }

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

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

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

633:     /*    fnet[2*lbus[i]]   += IDi; */
634:     row[0] = net_start + 2*lbus[i];
635:     col[0] = net_start + 2*lbus[i];  col[1] = net_start + 2*lbus[i]+1;
636:     val[0] = dIDi_dVr;               val[1] = dIDi_dVi;
637:     MatSetValues(J,1,row,2,col,val,ADD_VALUES);
638:     /*    fnet[2*lbus[i]+1] += IDr; */
639:     row[0] = net_start + 2*lbus[i]+1;
640:     col[0] = net_start + 2*lbus[i];  col[1] = net_start + 2*lbus[i]+1;
641:     val[0] = dIDr_dVr;               val[1] = dIDr_dVi;
642:     MatSetValues(J,1,row,2,col,val,ADD_VALUES);
643:   }
644:   VecRestoreArray(user->V0,&v0);

646:   VecRestoreArray(Xgen,&xgen);
647:   VecRestoreArray(Xnet,&xnet);

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

651:   MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
652:   MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
653:   return(0);
654: }

656: /*
657:    J = [I, 0
658:         dg_dx, dg_dy]
659: */
660: PetscErrorCode AlgJacobian(SNES snes,Vec X,Mat A,Mat B,void *ctx)
661: {
663:   Userctx        *user=(Userctx*)ctx;

666:   ResidualJacobian(snes,X,A,B,ctx);
667:   MatSetOption(A,MAT_KEEP_NONZERO_PATTERN,PETSC_TRUE);
668:   MatZeroRowsIS(A,user->is_diff,1.0,NULL,NULL);
669:   return(0);
670: }

672: /*
673:    J = [a*I-df_dx, -df_dy
674:         dg_dx, dg_dy]
675: */

677: PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,Userctx *user)
678: {
680:   SNES           snes;
681:   PetscScalar    atmp = (PetscScalar) a;
682:   PetscInt       i,row;

685:   user->t = t;

687:   TSGetSNES(ts,&snes);
688:   ResidualJacobian(snes,X,A,B,user);
689:   for (i=0;i < ngen;i++) {
690:     row = 9*i;
691:     MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
692:     row  = 9*i+1;
693:     MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
694:     row  = 9*i+2;
695:     MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
696:     row  = 9*i+3;
697:     MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
698:     row  = 9*i+6;
699:     MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
700:     row  = 9*i+7;
701:     MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
702:     row  = 9*i+8;
703:     MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
704:   }
705:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
706:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
707:   return(0);
708: }

710: static PetscErrorCode CostIntegrand(TS ts,PetscReal t,Vec U,Vec R,Userctx *user)
711: {
712:   PetscErrorCode    ierr;
713:   PetscScalar       *r;
714:   const PetscScalar *u;
715:   PetscInt          idx;
716:   Vec               Xgen,Xnet;
717:   PetscScalar       *xgen;
718:   PetscInt          i;

721:   DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
722:   DMCompositeScatter(user->dmpgrid,U,Xgen,Xnet);

724:   VecGetArray(Xgen,&xgen);

726:   VecGetArrayRead(U,&u);
727:   VecGetArray(R,&r);
728:   r[0] = 0.;

730:   idx = 0;
731:   for (i=0;i<ngen;i++) {
732:     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);
733:     idx  += 9;
734:   }
735:   VecRestoreArray(R,&r);
736:   VecRestoreArrayRead(U,&u);
737:   DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);
738:   return(0);
739: }

741: static PetscErrorCode MonitorUpdateQ(TS ts,PetscInt stepnum,PetscReal time,Vec X,void *ctx0)
742: {
744:   Vec            C,*Y;
745:   PetscInt       Nr;
746:   PetscReal      h,theta;
747:   Userctx        *ctx=(Userctx*)ctx0;

750:   theta = 0.5;
751:   TSGetStages(ts,&Nr,&Y);
752:   TSGetTimeStep(ts,&h);
753:   VecDuplicate(ctx->vec_q,&C);
754:   /* compute integrals */
755:   if (stepnum>0) {
756:     CostIntegrand(ts,time,X,C,ctx);
757:     VecAXPY(ctx->vec_q,h*theta,C);
758:     CostIntegrand(ts,time+h*theta,Y[0],C,ctx);
759:     VecAXPY(ctx->vec_q,h*(1-theta),C);
760:   }
761:   VecDestroy(&C);
762:   return(0);
763: }

765: int main(int argc,char **argv)
766: {
767:   Userctx            user;
768:   Vec                p;
769:   PetscScalar        *x_ptr;
770:   PetscErrorCode     ierr;
771:   PetscMPIInt        size;
772:   PetscInt           i;
773:   KSP                ksp;
774:   PC                 pc;
775:   PetscInt           *idx2;
776:   Tao                tao;
777:   Vec                lowerb,upperb;

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

784:   VecCreateSeq(PETSC_COMM_WORLD,1,&user.vec_q);

786:   user.neqs_gen   = 9*ngen; /* # eqs. for generator subsystem */
787:   user.neqs_net   = 2*nbus; /* # eqs. for network subsystem   */
788:   user.neqs_pgrid = user.neqs_gen + user.neqs_net;

790:   /* Create indices for differential and algebraic equations */
791:   PetscMalloc1(7*ngen,&idx2);
792:   for (i=0; i<ngen; i++) {
793:     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;
794:     idx2[7*i+4] = 9*i+6; idx2[7*i+5] = 9*i+7; idx2[7*i+6] = 9*i+8;
795:   }
796:   ISCreateGeneral(PETSC_COMM_WORLD,7*ngen,idx2,PETSC_COPY_VALUES,&user.is_diff);
797:   ISComplement(user.is_diff,0,user.neqs_pgrid,&user.is_alg);
798:   PetscFree(idx2);

800:   /* Set run time options */
801:   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Transient stability fault options","");
802:   {
803:     user.tfaulton  = 1.0;
804:     user.tfaultoff = 1.2;
805:     user.Rfault    = 0.0001;
806:     user.faultbus  = 8;
807:     PetscOptionsReal("-tfaulton","","",user.tfaulton,&user.tfaulton,NULL);
808:     PetscOptionsReal("-tfaultoff","","",user.tfaultoff,&user.tfaultoff,NULL);
809:     PetscOptionsInt("-faultbus","","",user.faultbus,&user.faultbus,NULL);
810:     user.t0        = 0.0;
811:     user.tmax      = 1.5;
812:     PetscOptionsReal("-t0","","",user.t0,&user.t0,NULL);
813:     PetscOptionsReal("-tmax","","",user.tmax,&user.tmax,NULL);
814:     user.freq_u    = 61.0;
815:     user.freq_l    = 59.0;
816:     user.pow       = 2;
817:     PetscOptionsReal("-frequ","","",user.freq_u,&user.freq_u,NULL);
818:     PetscOptionsReal("-freql","","",user.freq_l,&user.freq_l,NULL);
819:     PetscOptionsInt("-pow","","",user.pow,&user.pow,NULL);

821:   }
822:   PetscOptionsEnd();

824:   /* Create DMs for generator and network subsystems */
825:   DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,user.neqs_gen,1,1,NULL,&user.dmgen);
826:   DMSetOptionsPrefix(user.dmgen,"dmgen_");
827:   DMSetFromOptions(user.dmgen);
828:   DMSetUp(user.dmgen);
829:   DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,user.neqs_net,1,1,NULL,&user.dmnet);
830:   DMSetOptionsPrefix(user.dmnet,"dmnet_");
831:   DMSetFromOptions(user.dmnet);
832:   DMSetUp(user.dmnet);
833:   /* Create a composite DM packer and add the two DMs */
834:   DMCompositeCreate(PETSC_COMM_WORLD,&user.dmpgrid);
835:   DMSetOptionsPrefix(user.dmpgrid,"pgrid_");
836:   DMCompositeAddDM(user.dmpgrid,user.dmgen);
837:   DMCompositeAddDM(user.dmpgrid,user.dmnet);

839:   /* Create TAO solver and set desired solution method */
840:   TaoCreate(PETSC_COMM_WORLD,&tao);
841:   TaoSetType(tao,TAOBLMVM);
842:   /*
843:      Optimization starts
844:   */
845:   /* Set initial solution guess */
846:   VecCreateSeq(PETSC_COMM_WORLD,3,&p);
847:   VecGetArray(p,&x_ptr);
848:   x_ptr[0] = PG[0]; x_ptr[1] = PG[1]; x_ptr[2] = PG[2];
849:   VecRestoreArray(p,&x_ptr);

851:   TaoSetInitialVector(tao,p);
852:   /* Set routine for function and gradient evaluation */
853:   TaoSetObjectiveRoutine(tao,FormFunction,(void *)&user);
854:   TaoSetGradientRoutine(tao,TaoDefaultComputeGradient,(void *)&user);

856:   /* Set bounds for the optimization */
857:   VecDuplicate(p,&lowerb);
858:   VecDuplicate(p,&upperb);
859:   VecGetArray(lowerb,&x_ptr);
860:   x_ptr[0] = 0.5; x_ptr[1] = 0.5; x_ptr[2] = 0.5;
861:   VecRestoreArray(lowerb,&x_ptr);
862:   VecGetArray(upperb,&x_ptr);
863:   x_ptr[0] = 2.0; x_ptr[1] = 2.0; x_ptr[2] = 2.0;
864:   VecRestoreArray(upperb,&x_ptr);
865:   TaoSetVariableBounds(tao,lowerb,upperb);

867:   /* Check for any TAO command line options */
868:   TaoSetFromOptions(tao);
869:   TaoGetKSP(tao,&ksp);
870:   if (ksp) {
871:     KSPGetPC(ksp,&pc);
872:     PCSetType(pc,PCNONE);
873:   }

875:   /* SOLVE THE APPLICATION */
876:   TaoSolve(tao);

878:   VecView(p,PETSC_VIEWER_STDOUT_WORLD);
879:   /* Free TAO data structures */
880:   TaoDestroy(&tao);
881:   VecDestroy(&user.vec_q);
882:   VecDestroy(&lowerb);
883:   VecDestroy(&upperb);
884:   VecDestroy(&p);
885:   DMDestroy(&user.dmgen);
886:   DMDestroy(&user.dmnet);
887:   DMDestroy(&user.dmpgrid);
888:   ISDestroy(&user.is_diff);
889:   ISDestroy(&user.is_alg);
890:   PetscFinalize();
891:   return ierr;
892: }

894: /* ------------------------------------------------------------------ */
895: /*
896:    FormFunction - Evaluates the function and corresponding gradient.

898:    Input Parameters:
899:    tao - the Tao context
900:    X   - the input vector
901:    ptr - optional user-defined context, as set by TaoSetObjectiveAndGradientRoutine()

903:    Output Parameters:
904:    f   - the newly evaluated function
905: */
906: PetscErrorCode FormFunction(Tao tao,Vec P,PetscReal *f,void *ctx0)
907: {
908:   TS             ts;
909:   SNES           snes_alg;
911:   Userctx        *ctx = (Userctx*)ctx0;
912:   Vec            X;
913:   Mat            J;
914:   /* sensitivity context */
915:   PetscScalar    *x_ptr;
916:   PetscViewer    Xview,Ybusview;
917:   Vec            F_alg;
918:   Vec            Xdot;
919:   PetscInt       row_loc,col_loc;
920:   PetscScalar    val;

922:   VecGetArrayRead(P,(const PetscScalar**)&x_ptr);
923:   PG[0] = x_ptr[0];
924:   PG[1] = x_ptr[1];
925:   PG[2] = x_ptr[2];
926:   VecRestoreArrayRead(P,(const PetscScalar**)&x_ptr);

928:   ctx->stepnum = 0;

930:   VecZeroEntries(ctx->vec_q);

932:   /* Read initial voltage vector and Ybus */
933:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"X.bin",FILE_MODE_READ,&Xview);
934:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"Ybus.bin",FILE_MODE_READ,&Ybusview);

936:   VecCreate(PETSC_COMM_WORLD,&ctx->V0);
937:   VecSetSizes(ctx->V0,PETSC_DECIDE,ctx->neqs_net);
938:   VecLoad(ctx->V0,Xview);

940:   MatCreate(PETSC_COMM_WORLD,&ctx->Ybus);
941:   MatSetSizes(ctx->Ybus,PETSC_DECIDE,PETSC_DECIDE,ctx->neqs_net,ctx->neqs_net);
942:   MatSetType(ctx->Ybus,MATBAIJ);
943:   /*  MatSetBlockSize(ctx->Ybus,2); */
944:   MatLoad(ctx->Ybus,Ybusview);

946:   PetscViewerDestroy(&Xview);
947:   PetscViewerDestroy(&Ybusview);

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

951:   MatCreate(PETSC_COMM_WORLD,&J);
952:   MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,ctx->neqs_pgrid,ctx->neqs_pgrid);
953:   MatSetFromOptions(J);
954:   PreallocateJacobian(J,ctx);

956:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
957:      Create timestepping solver context
958:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
959:   TSCreate(PETSC_COMM_WORLD,&ts);
960:   TSSetProblemType(ts,TS_NONLINEAR);
961:   TSSetType(ts,TSCN);
962:   TSSetIFunction(ts,NULL,(TSIFunction) IFunction,ctx);
963:   TSSetIJacobian(ts,J,J,(TSIJacobian)IJacobian,ctx);
964:   TSSetApplicationContext(ts,ctx);

966:   TSMonitorSet(ts,MonitorUpdateQ,ctx,NULL);
967:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
968:      Set initial conditions
969:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
970:   SetInitialGuess(X,ctx);

972:   VecDuplicate(X,&F_alg);
973:   SNESCreate(PETSC_COMM_WORLD,&snes_alg);
974:   SNESSetFunction(snes_alg,F_alg,AlgFunction,ctx);
975:   MatZeroEntries(J);
976:   SNESSetJacobian(snes_alg,J,J,AlgJacobian,ctx);
977:   SNESSetOptionsPrefix(snes_alg,"alg_");
978:   SNESSetFromOptions(snes_alg);
979:   ctx->alg_flg = PETSC_TRUE;
980:   /* Solve the algebraic equations */
981:   SNESSolve(snes_alg,NULL,X);

983:   /* Just to set up the Jacobian structure */
984:   VecDuplicate(X,&Xdot);
985:   IJacobian(ts,0.0,X,Xdot,0.0,J,J,ctx);
986:   VecDestroy(&Xdot);

988:   ctx->stepnum++;

990:   TSSetTimeStep(ts,0.01);
991:   TSSetMaxTime(ts,ctx->tmax);
992:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);
993:   TSSetFromOptions(ts);

995:   /* Prefault period */
996:   ctx->alg_flg = PETSC_FALSE;
997:   TSSetTime(ts,0.0);
998:   TSSetMaxTime(ts,ctx->tfaulton);
999:   TSSolve(ts,X);

1001:   /* Create the nonlinear solver for solving the algebraic system */
1002:   /* Note that although the algebraic system needs to be solved only for
1003:      Idq and V, we reuse the entire system including xgen. The xgen
1004:      variables are held constant by setting their residuals to 0 and
1005:      putting a 1 on the Jacobian diagonal for xgen rows
1006:   */
1007:   MatZeroEntries(J);

1009:   /* Apply disturbance - resistive fault at ctx->faultbus */
1010:   /* This is done by adding shunt conductance to the diagonal location
1011:      in the Ybus matrix */
1012:   row_loc = 2*ctx->faultbus; col_loc = 2*ctx->faultbus+1; /* Location for G */
1013:   val     = 1/ctx->Rfault;
1014:   MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);
1015:   row_loc = 2*ctx->faultbus+1; col_loc = 2*ctx->faultbus; /* Location for G */
1016:   val     = 1/ctx->Rfault;
1017:   MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);

1019:   MatAssemblyBegin(ctx->Ybus,MAT_FINAL_ASSEMBLY);
1020:   MatAssemblyEnd(ctx->Ybus,MAT_FINAL_ASSEMBLY);

1022:   ctx->alg_flg = PETSC_TRUE;
1023:   /* Solve the algebraic equations */
1024:   SNESSolve(snes_alg,NULL,X);

1026:   ctx->stepnum++;

1028:   /* Disturbance period */
1029:   ctx->alg_flg = PETSC_FALSE;
1030:   TSSetTime(ts,ctx->tfaulton);
1031:   TSSetMaxTime(ts,ctx->tfaultoff);
1032:   TSSolve(ts,X);

1034:   /* Remove the fault */
1035:   row_loc = 2*ctx->faultbus; col_loc = 2*ctx->faultbus+1;
1036:   val     = -1/ctx->Rfault;
1037:   MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);
1038:   row_loc = 2*ctx->faultbus+1; col_loc = 2*ctx->faultbus;
1039:   val     = -1/ctx->Rfault;
1040:   MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);

1042:   MatAssemblyBegin(ctx->Ybus,MAT_FINAL_ASSEMBLY);
1043:   MatAssemblyEnd(ctx->Ybus,MAT_FINAL_ASSEMBLY);

1045:   MatZeroEntries(J);

1047:   ctx->alg_flg = PETSC_TRUE;

1049:   /* Solve the algebraic equations */
1050:   SNESSolve(snes_alg,NULL,X);

1052:   ctx->stepnum++;

1054:   /* Post-disturbance period */
1055:   ctx->alg_flg = PETSC_TRUE;
1056:   TSSetTime(ts,ctx->tfaultoff);
1057:   TSSetMaxTime(ts,ctx->tmax);
1058:   TSSolve(ts,X);

1060:   VecGetArray(ctx->vec_q,&x_ptr);
1061:   *f   = x_ptr[0];
1062:   VecRestoreArray(ctx->vec_q,&x_ptr);

1064:   MatDestroy(&ctx->Ybus);
1065:   VecDestroy(&ctx->V0);
1066:   SNESDestroy(&snes_alg);
1067:   VecDestroy(&F_alg);
1068:   MatDestroy(&J);
1069:   VecDestroy(&X);
1070:   TSDestroy(&ts);

1072:   return 0;
1073: }

1075: /*TEST

1077:   build:
1078:       requires: double !complex !defined(USE_64BIT_INDICES)

1080:    test:
1081:       args: -viewer_binary_skip_info -tao_monitor -tao_gttol .2
1082:       localrunfiles: petscoptions X.bin Ybus.bin

1084: TEST*/