Actual source code: ex9busopt_fd.c

petsc-3.7.5 2017-01-01
Report Typos and Errors
  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;


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

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

123: PetscErrorCode SetInitialGuess(Vec X,Userctx *user)
124: {
126:   Vec            Xgen,Xnet;
127:   PetscScalar    *xgen,*xnet;
128:   PetscInt       i,idx=0;
129:   PetscScalar    Vr,Vi,IGr,IGi,Vm,Vm2;
130:   PetscScalar    Eqp,Edp,delta;
131:   PetscScalar    Efd,RF,VR; /* Exciter variables */
132:   PetscScalar    Id,Iq;  /* Generator dq axis currents */
133:   PetscScalar    theta,Vd,Vq,SE;

136:   M[0] = 2*H[0]/w_s; M[1] = 2*H[1]/w_s; M[2] = 2*H[2]/w_s;
137:   D[0] = 0.1*M[0]; D[1] = 0.1*M[1]; D[2] = 0.1*M[2];

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

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

144:   /* Generator subsystem initialization */
145:   VecGetArray(Xgen,&xgen);
146:   VecGetArray(Xnet,&xnet);

148:   for (i=0; i < ngen; i++) {
149:     Vr  = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */
150:     Vi  = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */
151:     Vm  = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm;
152:     IGr = (Vr*PG[i] + Vi*QG[i])/Vm2;
153:     IGi = (Vi*PG[i] - Vr*QG[i])/Vm2;

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

157:     theta = PETSC_PI/2.0 - delta;

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

162:     Vd = Vr*PetscCosScalar(theta) - Vi*PetscSinScalar(theta);
163:     Vq = Vr*PetscSinScalar(theta) + Vi*PetscCosScalar(theta);

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

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

170:     /* The generator variables are ordered as [Eqp,Edp,delta,w,Id,Iq] */
171:     xgen[idx]   = Eqp;
172:     xgen[idx+1] = Edp;
173:     xgen[idx+2] = delta;
174:     xgen[idx+3] = w_s;

176:     idx = idx + 4;

178:     xgen[idx]   = Id;
179:     xgen[idx+1] = Iq;

181:     idx = idx + 2;

183:     /* Exciter */
184:     Efd = Eqp + (Xd[i] - Xdp[i])*Id;
185:     SE  = k1[i]*PetscExpScalar(k2[i]*Efd);
186:     VR  =  KE[i]*Efd + SE;
187:     RF  =  KF[i]*Efd/TF[i];

189:     xgen[idx]   = Efd;
190:     xgen[idx+1] = RF;
191:     xgen[idx+2] = VR;

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

195:     idx = idx + 3;
196:   }

198:   VecRestoreArray(Xgen,&xgen);
199:   VecRestoreArray(Xnet,&xnet);

201:   /* VecView(Xgen,0); */
202:   DMCompositeGather(user->dmpgrid,X,INSERT_VALUES,Xgen,Xnet);
203:   DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);
204:   return(0);
205: }

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

227:   VecZeroEntries(F);
228:   DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
229:   DMCompositeGetLocalVectors(user->dmpgrid,&Fgen,&Fnet);
230:   DMCompositeScatter(user->dmpgrid,X,Xgen,Xnet);
231:   DMCompositeScatter(user->dmpgrid,F,Fgen,Fnet);

233:   /* Network current balance residual IG + Y*V + IL = 0. Only YV is added here.
234:      The generator current injection, IG, and load current injection, ID are added later
235:   */
236:   /* Note that the values in Ybus are stored assuming the imaginary current balance
237:      equation is ordered first followed by real current balance equation for each bus.
238:      Thus imaginary current contribution goes in location 2*i, and
239:      real current contribution in 2*i+1
240:   */
241:   MatMult(user->Ybus,Xnet,Fnet);

243:   VecGetArray(Xgen,&xgen);
244:   VecGetArray(Xnet,&xnet);
245:   VecGetArray(Fgen,&fgen);
246:   VecGetArray(Fnet,&fnet);

248:   /* Generator subsystem */
249:   for (i=0; i < ngen; i++) {
250:     Eqp   = xgen[idx];
251:     Edp   = xgen[idx+1];
252:     delta = xgen[idx+2];
253:     w     = xgen[idx+3];
254:     Id    = xgen[idx+4];
255:     Iq    = xgen[idx+5];
256:     Efd   = xgen[idx+6];
257:     RF    = xgen[idx+7];
258:     VR    = xgen[idx+8];

260:     /* Generator differential equations */
261:     fgen[idx]   = (Eqp + (Xd[i] - Xdp[i])*Id - Efd)/Td0p[i];
262:     fgen[idx+1] = (Edp - (Xq[i] - Xqp[i])*Iq)/Tq0p[i];
263:     fgen[idx+2] = -w + w_s;
264:     fgen[idx+3] = (-TM[i] + Edp*Id + Eqp*Iq + (Xqp[i] - Xdp[i])*Id*Iq + D[i]*(w - w_s))/M[i];

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

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

273:     Zdq_inv[0] = Rs[i]/det;
274:     Zdq_inv[1] = Xqp[i]/det;
275:     Zdq_inv[2] = -Xdp[i]/det;
276:     Zdq_inv[3] = Rs[i]/det;

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

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

284:     fnet[2*gbus[i]]   -= IGi;
285:     fnet[2*gbus[i]+1] -= IGr;

287:     Vm = PetscSqrtScalar(Vd*Vd + Vq*Vq); Vm2 = Vm*Vm;

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

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

296:     idx = idx + 9;
297:   }

299:   VecGetArray(user->V0,&v0);
300:   for (i=0; i < nload; i++) {
301:     Vr  = xnet[2*lbus[i]]; /* Real part of load bus voltage */
302:     Vi  = xnet[2*lbus[i]+1]; /* Imaginary part of the load bus voltage */
303:     Vm  = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm;
304:     Vm0 = PetscSqrtScalar(v0[2*lbus[i]]*v0[2*lbus[i]] + v0[2*lbus[i]+1]*v0[2*lbus[i]+1]);
305:     PD  = QD = 0.0;
306:     for (k=0; k < ld_nsegsp[i]; k++) PD += ld_alphap[k]*PD0[i]*PetscPowScalar((Vm/Vm0),ld_betap[k]);
307:     for (k=0; k < ld_nsegsq[i]; k++) QD += ld_alphaq[k]*QD0[i]*PetscPowScalar((Vm/Vm0),ld_betaq[k]);

309:     /* Load currents */
310:     IDr = (PD*Vr + QD*Vi)/Vm2;
311:     IDi = (-QD*Vr + PD*Vi)/Vm2;

313:     fnet[2*lbus[i]]   += IDi;
314:     fnet[2*lbus[i]+1] += IDr;
315:   }
316:   VecRestoreArray(user->V0,&v0);

318:   VecRestoreArray(Xgen,&xgen);
319:   VecRestoreArray(Xnet,&xnet);
320:   VecRestoreArray(Fgen,&fgen);
321:   VecRestoreArray(Fnet,&fnet);

323:   DMCompositeGather(user->dmpgrid,F,INSERT_VALUES,Fgen,Fnet);
324:   DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);
325:   DMCompositeRestoreLocalVectors(user->dmpgrid,&Fgen,&Fnet);
326:   return(0);
327: }

329: /* \dot{x} - f(x,y)
330:      g(x,y) = 0
331:  */
334: PetscErrorCode IFunction(TS ts,PetscReal t, Vec X, Vec Xdot, Vec F, Userctx *user)
335: {
336:   PetscErrorCode    ierr;
337:   SNES              snes;
338:   PetscScalar       *f;
339:   const PetscScalar *xdot;
340:   PetscInt          i;

343:   user->t = t;

345:   TSGetSNES(ts,&snes);
346:   ResidualFunction(snes,X,F,user);
347:   VecGetArray(F,&f);
348:   VecGetArrayRead(Xdot,&xdot);
349:   for (i=0;i < ngen;i++) {
350:     f[9*i]   += xdot[9*i];
351:     f[9*i+1] += xdot[9*i+1];
352:     f[9*i+2] += xdot[9*i+2];
353:     f[9*i+3] += xdot[9*i+3];
354:     f[9*i+6] += xdot[9*i+6];
355:     f[9*i+7] += xdot[9*i+7];
356:     f[9*i+8] += xdot[9*i+8];
357:   }
358:   VecRestoreArray(F,&f);
359:   VecRestoreArrayRead(Xdot,&xdot);
360:   return(0);
361: }

363: /* This function is used for solving the algebraic system only during fault on and
364:    off times. It computes the entire F and then zeros out the part corresponding to
365:    differential equations
366:  F = [0;g(y)];
367: */
370: PetscErrorCode AlgFunction(SNES snes, Vec X, Vec F, void *ctx)
371: {
373:   Userctx        *user=(Userctx*)ctx;
374:   PetscScalar    *f;
375:   PetscInt       i;

378:   ResidualFunction(snes,X,F,user);
379:   VecGetArray(F,&f);
380:   for (i=0; i < ngen; i++) {
381:     f[9*i]   = 0;
382:     f[9*i+1] = 0;
383:     f[9*i+2] = 0;
384:     f[9*i+3] = 0;
385:     f[9*i+6] = 0;
386:     f[9*i+7] = 0;
387:     f[9*i+8] = 0;
388:   }
389:   VecRestoreArray(F,&f);
390:   return(0);
391: }

395: PetscErrorCode PreallocateJacobian(Mat J, Userctx *user)
396: {
398:   PetscInt       *d_nnz;
399:   PetscInt       i,idx=0,start=0;
400:   PetscInt       ncols;

403:   PetscMalloc1(user->neqs_pgrid,&d_nnz);
404:   for (i=0; i<user->neqs_pgrid; i++) d_nnz[i] = 0;
405:   /* Generator subsystem */
406:   for (i=0; i < ngen; i++) {

408:     d_nnz[idx]   += 3;
409:     d_nnz[idx+1] += 2;
410:     d_nnz[idx+2] += 2;
411:     d_nnz[idx+3] += 5;
412:     d_nnz[idx+4] += 6;
413:     d_nnz[idx+5] += 6;

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

418:     d_nnz[idx+6] += 2;
419:     d_nnz[idx+7] += 2;
420:     d_nnz[idx+8] += 5;

422:     idx = idx + 9;
423:   }

425:   start = user->neqs_gen;

427:   for (i=0; i < nbus; i++) {
428:     MatGetRow(user->Ybus,2*i,&ncols,NULL,NULL);
429:     d_nnz[start+2*i]   += ncols;
430:     d_nnz[start+2*i+1] += ncols;
431:     MatRestoreRow(user->Ybus,2*i,&ncols,NULL,NULL);
432:   }

434:   MatSeqAIJSetPreallocation(J,0,d_nnz);

436:   PetscFree(d_nnz);
437:   return(0);
438: }

440: /*
441:    J = [-df_dx, -df_dy
442:         dg_dx, dg_dy]
443: */
446: PetscErrorCode ResidualJacobian(SNES snes,Vec X,Mat J,Mat B,void *ctx)
447: {
448:   PetscErrorCode    ierr;
449:   Userctx           *user=(Userctx*)ctx;
450:   Vec               Xgen,Xnet;
451:   PetscScalar       *xgen,*xnet;
452:   PetscInt          i,idx=0;
453:   PetscScalar       Vr,Vi,Vm,Vm2;
454:   PetscScalar       Eqp,Edp,delta; /* Generator variables */
455:   PetscScalar       Efd; /* Exciter variables */
456:   PetscScalar       Id,Iq;  /* Generator dq axis currents */
457:   PetscScalar       Vd,Vq;
458:   PetscScalar       val[10];
459:   PetscInt          row[2],col[10];
460:   PetscInt          net_start=user->neqs_gen;
461:   PetscScalar       Zdq_inv[4],det;
462:   PetscScalar       dVd_dVr,dVd_dVi,dVq_dVr,dVq_dVi,dVd_ddelta,dVq_ddelta;
463:   PetscScalar       dIGr_ddelta,dIGi_ddelta,dIGr_dId,dIGr_dIq,dIGi_dId,dIGi_dIq;
464:   PetscScalar       dSE_dEfd;
465:   PetscScalar       dVm_dVd,dVm_dVq,dVm_dVr,dVm_dVi;
466:   PetscInt          ncols;
467:   const PetscInt    *cols;
468:   const PetscScalar *yvals;
469:   PetscInt          k;
470:   PetscScalar       PD,QD,Vm0,*v0,Vm4;
471:   PetscScalar       dPD_dVr,dPD_dVi,dQD_dVr,dQD_dVi;
472:   PetscScalar       dIDr_dVr,dIDr_dVi,dIDi_dVr,dIDi_dVi;

475:   MatZeroEntries(B);
476:   DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
477:   DMCompositeScatter(user->dmpgrid,X,Xgen,Xnet);

479:   VecGetArray(Xgen,&xgen);
480:   VecGetArray(Xnet,&xnet);

482:   /* Generator subsystem */
483:   for (i=0; i < ngen; i++) {
484:     Eqp   = xgen[idx];
485:     Edp   = xgen[idx+1];
486:     delta = xgen[idx+2];
487:     Id    = xgen[idx+4];
488:     Iq    = xgen[idx+5];
489:     Efd   = xgen[idx+6];

491:     /*    fgen[idx]   = (Eqp + (Xd[i] - Xdp[i])*Id - Efd)/Td0p[i]; */
492:     row[0] = idx;
493:     col[0] = idx;           col[1] = idx+4;          col[2] = idx+6;
494:     val[0] = 1/ Td0p[i]; val[1] = (Xd[i] - Xdp[i])/ Td0p[i]; val[2] = -1/Td0p[i];

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

498:     /*    fgen[idx+1] = (Edp - (Xq[i] - Xqp[i])*Iq)/Tq0p[i]; */
499:     row[0] = idx + 1;
500:     col[0] = idx + 1;       col[1] = idx+5;
501:     val[0] = 1/Tq0p[i]; val[1] = -(Xq[i] - Xqp[i])/Tq0p[i];
502:     MatSetValues(J,1,row,2,col,val,INSERT_VALUES);

504:     /*    fgen[idx+2] = - w + w_s; */
505:     row[0] = idx + 2;
506:     col[0] = idx + 2; col[1] = idx + 3;
507:     val[0] = 0;       val[1] = -1;
508:     MatSetValues(J,1,row,2,col,val,INSERT_VALUES);

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

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

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

522:     Zdq_inv[0] = Rs[i]/det;
523:     Zdq_inv[1] = Xqp[i]/det;
524:     Zdq_inv[2] = -Xdp[i]/det;
525:     Zdq_inv[3] = Rs[i]/det;

527:     dVd_dVr    = PetscSinScalar(delta); dVd_dVi = -PetscCosScalar(delta);
528:     dVq_dVr    = PetscCosScalar(delta); dVq_dVi = PetscSinScalar(delta);
529:     dVd_ddelta = Vr*PetscCosScalar(delta) + Vi*PetscSinScalar(delta);
530:     dVq_ddelta = -Vr*PetscSinScalar(delta) + Vi*PetscCosScalar(delta);

532:     /*    fgen[idx+4] = Zdq_inv[0]*(-Edp + Vd) + Zdq_inv[1]*(-Eqp + Vq) + Id; */
533:     row[0] = idx+4;
534:     col[0] = idx;         col[1] = idx+1;        col[2] = idx + 2;
535:     val[0] = -Zdq_inv[1]; val[1] = -Zdq_inv[0];  val[2] = Zdq_inv[0]*dVd_ddelta + Zdq_inv[1]*dVq_ddelta;
536:     col[3] = idx + 4; col[4] = net_start+2*gbus[i];                     col[5] = net_start + 2*gbus[i]+1;
537:     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;
538:     MatSetValues(J,1,row,6,col,val,INSERT_VALUES);

540:     /*  fgen[idx+5] = Zdq_inv[2]*(-Edp + Vd) + Zdq_inv[3]*(-Eqp + Vq) + Iq; */
541:     row[0] = idx+5;
542:     col[0] = idx;         col[1] = idx+1;        col[2] = idx + 2;
543:     val[0] = -Zdq_inv[3]; val[1] = -Zdq_inv[2];  val[2] = Zdq_inv[2]*dVd_ddelta + Zdq_inv[3]*dVq_ddelta;
544:     col[3] = idx + 5; col[4] = net_start+2*gbus[i];                     col[5] = net_start + 2*gbus[i]+1;
545:     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;
546:     MatSetValues(J,1,row,6,col,val,INSERT_VALUES);

548:     dIGr_ddelta = Id*PetscCosScalar(delta) - Iq*PetscSinScalar(delta);
549:     dIGi_ddelta = Id*PetscSinScalar(delta) + Iq*PetscCosScalar(delta);
550:     dIGr_dId    = PetscSinScalar(delta);  dIGr_dIq = PetscCosScalar(delta);
551:     dIGi_dId    = -PetscCosScalar(delta); dIGi_dIq = PetscSinScalar(delta);

553:     /* fnet[2*gbus[i]]   -= IGi; */
554:     row[0] = net_start + 2*gbus[i];
555:     col[0] = idx+2;        col[1] = idx + 4;   col[2] = idx + 5;
556:     val[0] = -dIGi_ddelta; val[1] = -dIGi_dId; val[2] = -dIGi_dIq;
557:     MatSetValues(J,1,row,3,col,val,INSERT_VALUES);

559:     /* fnet[2*gbus[i]+1]   -= IGr; */
560:     row[0] = net_start + 2*gbus[i]+1;
561:     col[0] = idx+2;        col[1] = idx + 4;   col[2] = idx + 5;
562:     val[0] = -dIGr_ddelta; val[1] = -dIGr_dId; val[2] = -dIGr_dIq;
563:     MatSetValues(J,1,row,3,col,val,INSERT_VALUES);

565:     Vm = PetscSqrtScalar(Vd*Vd + Vq*Vq); Vm2 = Vm*Vm;

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

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

572:     row[0] = idx + 6;
573:     col[0] = idx + 6;                     col[1] = idx + 8;
574:     val[0] = (KE[i] + dSE_dEfd)/TE[i];  val[1] = -1/TE[i];
575:     MatSetValues(J,1,row,2,col,val,INSERT_VALUES);

577:     /* Exciter differential equations */

579:     /*    fgen[idx+7] = (RF - KF[i]*Efd/TF[i])/TF[i]; */
580:     row[0] = idx + 7;
581:     col[0] = idx + 6;       col[1] = idx + 7;
582:     val[0] = (-KF[i]/TF[i])/TF[i];  val[1] = 1/TF[i];
583:     MatSetValues(J,1,row,2,col,val,INSERT_VALUES);

585:     /*    fgen[idx+8] = (VR - KA[i]*RF + KA[i]*KF[i]*Efd/TF[i] - KA[i]*(Vref[i] - Vm))/TA[i]; */
586:     /* Vm = (Vd^2 + Vq^2)^0.5; */
587:     dVm_dVd    = Vd/Vm; dVm_dVq = Vq/Vm;
588:     dVm_dVr    = dVm_dVd*dVd_dVr + dVm_dVq*dVq_dVr;
589:     dVm_dVi    = dVm_dVd*dVd_dVi + dVm_dVq*dVq_dVi;
590:     row[0]     = idx + 8;
591:     col[0]     = idx + 6;           col[1] = idx + 7; col[2] = idx + 8;
592:     val[0]     = (KA[i]*KF[i]/TF[i])/TA[i]; val[1] = -KA[i]/TA[i];  val[2] = 1/TA[i];
593:     col[3]     = net_start + 2*gbus[i]; col[4] = net_start + 2*gbus[i]+1;
594:     val[3]     = KA[i]*dVm_dVr/TA[i];         val[4] = KA[i]*dVm_dVi/TA[i];
595:     MatSetValues(J,1,row,5,col,val,INSERT_VALUES);
596:     idx        = idx + 9;
597:   }


600:   for (i=0; i<nbus; i++) {
601:     MatGetRow(user->Ybus,2*i,&ncols,&cols,&yvals);
602:     row[0] = net_start + 2*i;
603:     for (k=0; k<ncols; k++) {
604:       col[k] = net_start + cols[k];
605:       val[k] = yvals[k];
606:     }
607:     MatSetValues(J,1,row,ncols,col,val,INSERT_VALUES);
608:     MatRestoreRow(user->Ybus,2*i,&ncols,&cols,&yvals);

610:     MatGetRow(user->Ybus,2*i+1,&ncols,&cols,&yvals);
611:     row[0] = net_start + 2*i+1;
612:     for (k=0; k<ncols; k++) {
613:       col[k] = net_start + cols[k];
614:       val[k] = yvals[k];
615:     }
616:     MatSetValues(J,1,row,ncols,col,val,INSERT_VALUES);
617:     MatRestoreRow(user->Ybus,2*i+1,&ncols,&cols,&yvals);
618:   }

620:   MatAssemblyBegin(J,MAT_FLUSH_ASSEMBLY);
621:   MatAssemblyEnd(J,MAT_FLUSH_ASSEMBLY);

623:   VecGetArray(user->V0,&v0);
624:   for (i=0; i < nload; i++) {
625:     Vr      = xnet[2*lbus[i]]; /* Real part of load bus voltage */
626:     Vi      = xnet[2*lbus[i]+1]; /* Imaginary part of the load bus voltage */
627:     Vm      = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm; Vm4 = Vm2*Vm2;
628:     Vm0     = PetscSqrtScalar(v0[2*lbus[i]]*v0[2*lbus[i]] + v0[2*lbus[i]+1]*v0[2*lbus[i]+1]);
629:     PD      = QD = 0.0;
630:     dPD_dVr = dPD_dVi = dQD_dVr = dQD_dVi = 0.0;
631:     for (k=0; k < ld_nsegsp[i]; k++) {
632:       PD      += ld_alphap[k]*PD0[i]*PetscPowScalar((Vm/Vm0),ld_betap[k]);
633:       dPD_dVr += ld_alphap[k]*ld_betap[k]*PD0[i]*PetscPowScalar((1/Vm0),ld_betap[k])*Vr*PetscPowScalar(Vm,(ld_betap[k]-2));
634:       dPD_dVi += ld_alphap[k]*ld_betap[k]*PD0[i]*PetscPowScalar((1/Vm0),ld_betap[k])*Vi*PetscPowScalar(Vm,(ld_betap[k]-2));
635:     }
636:     for (k=0; k < ld_nsegsq[i]; k++) {
637:       QD      += ld_alphaq[k]*QD0[i]*PetscPowScalar((Vm/Vm0),ld_betaq[k]);
638:       dQD_dVr += ld_alphaq[k]*ld_betaq[k]*QD0[i]*PetscPowScalar((1/Vm0),ld_betaq[k])*Vr*PetscPowScalar(Vm,(ld_betaq[k]-2));
639:       dQD_dVi += ld_alphaq[k]*ld_betaq[k]*QD0[i]*PetscPowScalar((1/Vm0),ld_betaq[k])*Vi*PetscPowScalar(Vm,(ld_betaq[k]-2));
640:     }

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

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

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


652:     /*    fnet[2*lbus[i]]   += IDi; */
653:     row[0] = net_start + 2*lbus[i];
654:     col[0] = net_start + 2*lbus[i];  col[1] = net_start + 2*lbus[i]+1;
655:     val[0] = dIDi_dVr;               val[1] = dIDi_dVi;
656:     MatSetValues(J,1,row,2,col,val,ADD_VALUES);
657:     /*    fnet[2*lbus[i]+1] += IDr; */
658:     row[0] = net_start + 2*lbus[i]+1;
659:     col[0] = net_start + 2*lbus[i];  col[1] = net_start + 2*lbus[i]+1;
660:     val[0] = dIDr_dVr;               val[1] = dIDr_dVi;
661:     MatSetValues(J,1,row,2,col,val,ADD_VALUES);
662:   }
663:   VecRestoreArray(user->V0,&v0);

665:   VecRestoreArray(Xgen,&xgen);
666:   VecRestoreArray(Xnet,&xnet);

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

670:   MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
671:   MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
672:   return(0);
673: }

675: /*
676:    J = [I, 0
677:         dg_dx, dg_dy]
678: */
681: PetscErrorCode AlgJacobian(SNES snes,Vec X,Mat A,Mat B,void *ctx)
682: {
684:   Userctx        *user=(Userctx*)ctx;

687:   ResidualJacobian(snes,X,A,B,ctx);
688:   MatSetOption(A,MAT_KEEP_NONZERO_PATTERN,PETSC_TRUE);
689:   MatZeroRowsIS(A,user->is_diff,1.0,NULL,NULL);
690:   return(0);
691: }

693: /*
694:    J = [a*I-df_dx, -df_dy
695:         dg_dx, dg_dy]
696: */

700: PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,Userctx *user)
701: {
703:   SNES           snes;
704:   PetscScalar    atmp = (PetscScalar) a;
705:   PetscInt       i,row;

708:   user->t = t;

710:   TSGetSNES(ts,&snes);
711:   ResidualJacobian(snes,X,A,B,user);
712:   for (i=0;i < ngen;i++) {
713:     row = 9*i;
714:     MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
715:     row  = 9*i+1;
716:     MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
717:     row  = 9*i+2;
718:     MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
719:     row  = 9*i+3;
720:     MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
721:     row  = 9*i+6;
722:     MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
723:     row  = 9*i+7;
724:     MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
725:     row  = 9*i+8;
726:     MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES);
727:   }
728:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
729:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
730:   return(0);
731: }

735: static PetscErrorCode CostIntegrand(TS ts,PetscReal t,Vec U,Vec R,Userctx *user)
736: {
737:   PetscErrorCode    ierr;
738:   PetscScalar       *r;
739:   const PetscScalar *u;
740:   PetscInt          idx;
741:   Vec               Xgen,Xnet;
742:   PetscScalar       *xgen;
743:   PetscInt          i;

746:   DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
747:   DMCompositeScatter(user->dmpgrid,U,Xgen,Xnet);

749:   VecGetArray(Xgen,&xgen);

751:   VecGetArrayRead(U,&u);
752:   VecGetArray(R,&r);
753:   r[0] = 0.;

755:   idx = 0;
756:   for (i=0;i<ngen;i++) {
757:     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);
758:     idx  += 9;
759:   }
760:   VecRestoreArray(R,&r);
761:   VecRestoreArrayRead(U,&u);
762:   DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet);
763:   return(0);
764: }

768: static PetscErrorCode MonitorUpdateQ(TS ts,PetscInt stepnum,PetscReal time,Vec X,void *ctx0)
769: {
771:   Vec            C,*Y;
772:   PetscInt       Nr;
773:   PetscReal      h,theta;
774:   Userctx        *ctx=(Userctx*)ctx0;

777:   theta = 0.5;
778:   TSGetStages(ts,&Nr,&Y);
779:   TSGetTimeStep(ts,&h);
780:   VecDuplicate(ctx->vec_q,&C);
781:   /* compute integrals */
782:   if (stepnum>0) {
783:     CostIntegrand(ts,time,X,C,ctx);
784:     VecAXPY(ctx->vec_q,h*theta,C);
785:     CostIntegrand(ts,time+h*theta,Y[0],C,ctx);
786:     VecAXPY(ctx->vec_q,h*(1-theta),C);
787:   }
788:   VecDestroy(&C);
789:   return(0);
790: }

794: int main(int argc,char **argv)
795: {
796:   Userctx            user;
797:   Vec                p;
798:   PetscScalar        *x_ptr;
799:   PetscErrorCode     ierr;
800:   PetscMPIInt        size;
801:   PetscInt           i;
802:   KSP                ksp;
803:   PC                 pc;
804:   PetscInt           *idx2;
805:   Tao                tao;
806:   Vec                lowerb,upperb;

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

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

815:   user.neqs_gen   = 9*ngen; /* # eqs. for generator subsystem */
816:   user.neqs_net   = 2*nbus; /* # eqs. for network subsystem   */
817:   user.neqs_pgrid = user.neqs_gen + user.neqs_net;

819:   /* Create indices for differential and algebraic equations */
820:   PetscMalloc1(7*ngen,&idx2);
821:   for (i=0; i<ngen; i++) {
822:     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;
823:     idx2[7*i+4] = 9*i+6; idx2[7*i+5] = 9*i+7; idx2[7*i+6] = 9*i+8;
824:   }
825:   ISCreateGeneral(PETSC_COMM_WORLD,7*ngen,idx2,PETSC_COPY_VALUES,&user.is_diff);
826:   ISComplement(user.is_diff,0,user.neqs_pgrid,&user.is_alg);
827:   PetscFree(idx2);

829:   /* Set run time options */
830:   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Transient stability fault options","");
831:   {
832:     user.tfaulton  = 1.0;
833:     user.tfaultoff = 1.2;
834:     user.Rfault    = 0.0001;
835:     user.faultbus  = 8;
836:     PetscOptionsReal("-tfaulton","","",user.tfaulton,&user.tfaulton,NULL);
837:     PetscOptionsReal("-tfaultoff","","",user.tfaultoff,&user.tfaultoff,NULL);
838:     PetscOptionsInt("-faultbus","","",user.faultbus,&user.faultbus,NULL);
839:     user.t0        = 0.0;
840:     user.tmax      = 1.5;
841:     PetscOptionsReal("-t0","","",user.t0,&user.t0,NULL);
842:     PetscOptionsReal("-tmax","","",user.tmax,&user.tmax,NULL);
843:     user.freq_u    = 61.0;
844:     user.freq_l    = 59.0;
845:     user.pow       = 2;
846:     PetscOptionsReal("-frequ","","",user.freq_u,&user.freq_u,NULL);
847:     PetscOptionsReal("-freql","","",user.freq_l,&user.freq_l,NULL);
848:     PetscOptionsInt("-pow","","",user.pow,&user.pow,NULL);

850:   }
851:   PetscOptionsEnd();

853:   /* Create DMs for generator and network subsystems */
854:   DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,user.neqs_gen,1,1,NULL,&user.dmgen);
855:   DMSetOptionsPrefix(user.dmgen,"dmgen_");
856:   DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,user.neqs_net,1,1,NULL,&user.dmnet);
857:   DMSetOptionsPrefix(user.dmnet,"dmnet_");
858:   /* Create a composite DM packer and add the two DMs */
859:   DMCompositeCreate(PETSC_COMM_WORLD,&user.dmpgrid);
860:   DMSetOptionsPrefix(user.dmpgrid,"pgrid_");
861:   DMCompositeAddDM(user.dmpgrid,user.dmgen);
862:   DMCompositeAddDM(user.dmpgrid,user.dmnet);

864:   /* Create TAO solver and set desired solution method */
865:   TaoCreate(PETSC_COMM_WORLD,&tao);
866:   TaoSetType(tao,TAOBLMVM);
867:   /*
868:      Optimization starts
869:   */
870:   /* Set initial solution guess */
871:   VecCreateSeq(PETSC_COMM_WORLD,3,&p);
872:   VecGetArray(p,&x_ptr);
873:   x_ptr[0] = PG[0]; x_ptr[1] = PG[1]; x_ptr[2] = PG[2];
874:   VecRestoreArray(p,&x_ptr);

876:   TaoSetInitialVector(tao,p);
877:   /* Set routine for function and gradient evaluation */
878:   TaoSetObjectiveRoutine(tao,FormFunction,(void *)&user);
879:   TaoSetGradientRoutine(tao,TaoDefaultComputeGradient,(void *)&user);

881:   /* Set bounds for the optimization */
882:   VecDuplicate(p,&lowerb);
883:   VecDuplicate(p,&upperb);
884:   VecGetArray(lowerb,&x_ptr);
885:   x_ptr[0] = 0.5; x_ptr[1] = 0.5; x_ptr[2] = 0.5;
886:   VecRestoreArray(lowerb,&x_ptr);
887:   VecGetArray(upperb,&x_ptr);
888:   x_ptr[0] = 2.0; x_ptr[1] = 2.0; x_ptr[2] = 2.0;
889:   VecRestoreArray(upperb,&x_ptr);
890:   TaoSetVariableBounds(tao,lowerb,upperb);

892:   /* Check for any TAO command line options */
893:   TaoSetFromOptions(tao);
894:   TaoGetKSP(tao,&ksp);
895:   if (ksp) {
896:     KSPGetPC(ksp,&pc);
897:     PCSetType(pc,PCNONE);
898:   }

900:   /* SOLVE THE APPLICATION */
901:   TaoSolve(tao);

903:   VecView(p,PETSC_VIEWER_STDOUT_WORLD);
904:   /* Free TAO data structures */
905:   TaoDestroy(&tao);
906:   VecDestroy(&user.vec_q);
907:   VecDestroy(&lowerb);
908:   VecDestroy(&upperb);
909:   VecDestroy(&p);
910:   DMDestroy(&user.dmgen);
911:   DMDestroy(&user.dmnet);
912:   DMDestroy(&user.dmpgrid);
913:   ISDestroy(&user.is_diff);
914:   ISDestroy(&user.is_alg);
915:   PetscFinalize();
916:   return(0);
917: }

919: /* ------------------------------------------------------------------ */
922: /*
923:    FormFunction - Evaluates the function and corresponding gradient.

925:    Input Parameters:
926:    tao - the Tao context
927:    X   - the input vector
928:    ptr - optional user-defined context, as set by TaoSetObjectiveAndGradientRoutine()

930:    Output Parameters:
931:    f   - the newly evaluated function
932: */
933: PetscErrorCode FormFunction(Tao tao,Vec P,PetscReal *f,void *ctx0)
934: {
935:   TS             ts;
936:   SNES           snes_alg;
938:   Userctx        *ctx = (Userctx*)ctx0;
939:   Vec            X;
940:   Mat            J;
941:   /* sensitivity context */
942:   PetscScalar    *x_ptr;
943:   PetscViewer    Xview,Ybusview;
944:   Vec            F_alg;
945:   Vec            Xdot;
946:   PetscInt       row_loc,col_loc;
947:   PetscScalar    val;

949:   VecGetArray(P,&x_ptr);
950:   PG[0] = x_ptr[0];
951:   PG[1] = x_ptr[1];
952:   PG[2] = x_ptr[2];
953:   VecRestoreArray(P,&x_ptr);

955:   ctx->stepnum = 0;

957:   VecZeroEntries(ctx->vec_q);

959:   /* Read initial voltage vector and Ybus */
960:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"X.bin",FILE_MODE_READ,&Xview);
961:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"Ybus.bin",FILE_MODE_READ,&Ybusview);

963:   VecCreate(PETSC_COMM_WORLD,&ctx->V0);
964:   VecSetSizes(ctx->V0,PETSC_DECIDE,ctx->neqs_net);
965:   VecLoad(ctx->V0,Xview);

967:   MatCreate(PETSC_COMM_WORLD,&ctx->Ybus);
968:   MatSetSizes(ctx->Ybus,PETSC_DECIDE,PETSC_DECIDE,ctx->neqs_net,ctx->neqs_net);
969:   MatSetType(ctx->Ybus,MATBAIJ);
970:   /*  MatSetBlockSize(ctx->Ybus,2); */
971:   MatLoad(ctx->Ybus,Ybusview);

973:   PetscViewerDestroy(&Xview);
974:   PetscViewerDestroy(&Ybusview);

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

978:   MatCreate(PETSC_COMM_WORLD,&J);
979:   MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,ctx->neqs_pgrid,ctx->neqs_pgrid);
980:   MatSetFromOptions(J);
981:   PreallocateJacobian(J,ctx);

983:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
984:      Create timestepping solver context
985:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
986:   TSCreate(PETSC_COMM_WORLD,&ts);
987:   TSSetProblemType(ts,TS_NONLINEAR);
988:   TSSetType(ts,TSCN);
989:   TSSetIFunction(ts,NULL,(TSIFunction) IFunction,ctx);
990:   TSSetIJacobian(ts,J,J,(TSIJacobian)IJacobian,ctx);
991:   TSSetApplicationContext(ts,ctx);

993:   TSMonitorSet(ts,MonitorUpdateQ,ctx,NULL);
994:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
995:      Set initial conditions
996:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
997:   SetInitialGuess(X,ctx);

999:   VecDuplicate(X,&F_alg);
1000:   SNESCreate(PETSC_COMM_WORLD,&snes_alg);
1001:   SNESSetFunction(snes_alg,F_alg,AlgFunction,ctx);
1002:   MatZeroEntries(J);
1003:   SNESSetJacobian(snes_alg,J,J,AlgJacobian,ctx);
1004:   SNESSetOptionsPrefix(snes_alg,"alg_");
1005:   SNESSetFromOptions(snes_alg);
1006:   ctx->alg_flg = PETSC_TRUE;
1007:   /* Solve the algebraic equations */
1008:   SNESSolve(snes_alg,NULL,X);

1010:   /* Just to set up the Jacobian structure */
1011:   VecDuplicate(X,&Xdot);
1012:   IJacobian(ts,0.0,X,Xdot,0.0,J,J,ctx);
1013:   VecDestroy(&Xdot);

1015:   ctx->stepnum++;

1017:   TSSetDuration(ts,1000,ctx->tfaulton);
1018:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
1019:   TSSetInitialTimeStep(ts,0.0,0.01);
1020:   TSSetFromOptions(ts);
1021:   /* TSSetPostStep(ts,SaveSolution); */

1023:   ctx->alg_flg = PETSC_FALSE;
1024:   /* Prefault period */
1025:   TSSolve(ts,X);

1027:   /* Create the nonlinear solver for solving the algebraic system */
1028:   /* Note that although the algebraic system needs to be solved only for
1029:      Idq and V, we reuse the entire system including xgen. The xgen
1030:      variables are held constant by setting their residuals to 0 and
1031:      putting a 1 on the Jacobian diagonal for xgen rows
1032:   */
1033:   MatZeroEntries(J);

1035:   /* Apply disturbance - resistive fault at ctx->faultbus */
1036:   /* This is done by adding shunt conductance to the diagonal location
1037:      in the Ybus matrix */
1038:   row_loc = 2*ctx->faultbus; col_loc = 2*ctx->faultbus+1; /* Location for G */
1039:   val     = 1/ctx->Rfault;
1040:   MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);
1041:   row_loc = 2*ctx->faultbus+1; col_loc = 2*ctx->faultbus; /* Location for G */
1042:   val     = 1/ctx->Rfault;
1043:   MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);

1045:   MatAssemblyBegin(ctx->Ybus,MAT_FINAL_ASSEMBLY);
1046:   MatAssemblyEnd(ctx->Ybus,MAT_FINAL_ASSEMBLY);

1048:   ctx->alg_flg = PETSC_TRUE;
1049:   /* Solve the algebraic equations */
1050:   SNESSolve(snes_alg,NULL,X);

1052:   ctx->stepnum++;

1054:   /* Disturbance period */
1055:   TSSetDuration(ts,1000,ctx->tfaultoff);
1056:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
1057:   TSSetInitialTimeStep(ts,ctx->tfaulton,.01);

1059:   ctx->alg_flg = PETSC_FALSE;

1061:   TSSolve(ts,X);

1063:   /* Remove the fault */
1064:   row_loc = 2*ctx->faultbus; col_loc = 2*ctx->faultbus+1;
1065:   val     = -1/ctx->Rfault;
1066:   MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);
1067:   row_loc = 2*ctx->faultbus+1; col_loc = 2*ctx->faultbus;
1068:   val     = -1/ctx->Rfault;
1069:   MatSetValues(ctx->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES);

1071:   MatAssemblyBegin(ctx->Ybus,MAT_FINAL_ASSEMBLY);
1072:   MatAssemblyEnd(ctx->Ybus,MAT_FINAL_ASSEMBLY);

1074:   MatZeroEntries(J);

1076:   ctx->alg_flg = PETSC_TRUE;

1078:   /* Solve the algebraic equations */
1079:   SNESSolve(snes_alg,NULL,X);

1081:   ctx->stepnum++;

1083:   /* Post-disturbance period */
1084:   TSSetDuration(ts,1000,ctx->tmax);
1085:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
1086:   TSSetInitialTimeStep(ts,ctx->tfaultoff,.01);

1088:   ctx->alg_flg = PETSC_TRUE;

1090:   TSSolve(ts,X);
1091:   VecGetArray(ctx->vec_q,&x_ptr);
1092:   *f   = x_ptr[0];
1093:   VecRestoreArray(ctx->vec_q,&x_ptr);

1095:   MatDestroy(&ctx->Ybus);
1096:   VecDestroy(&ctx->V0);
1097:   SNESDestroy(&snes_alg);
1098:   VecDestroy(&F_alg);
1099:   MatDestroy(&J);
1100:   VecDestroy(&X);
1101:   TSDestroy(&ts);

1103:   return 0;
1104: }