Actual source code: ex9busoptfd.c

petsc-3.13.4 2020-08-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 */
100: PetscErrorCode dq2ri(PetscScalar Fd,PetscScalar Fq,PetscScalar delta,PetscScalar *Fr, PetscScalar *Fi)
101: {
103:   *Fr =  Fd*PetscSinScalar(delta) + Fq*PetscCosScalar(delta);
104:   *Fi = -Fd*PetscCosScalar(delta) + Fq*PetscSinScalar(delta);
105:   return(0);
106: }

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

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

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

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

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

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

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

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

151:     theta = PETSC_PI/2.0 - delta;

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

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

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

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

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

170:     idx = idx + 4;

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

175:     idx = idx + 2;

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

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

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

189:     idx = idx + 3;
190:   }

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

288:     idx = idx + 9;
289:   }

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

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

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

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

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

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

333:   user->t = t;

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

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

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

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

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

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

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

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

408:     idx = idx + 9;
409:   }

411:   start = user->neqs_gen;

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

420:   MatSeqAIJSetPreallocation(J,0,d_nnz);

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

561:     /* Exciter differential equations */

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

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


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

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

604:   MatAssemblyBegin(J,MAT_FLUSH_ASSEMBLY);
605:   MatAssemblyEnd(J,MAT_FLUSH_ASSEMBLY);

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

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

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

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


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

649:   VecRestoreArray(Xgen,&xgen);
650:   VecRestoreArray(Xnet,&xnet);

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

654:   MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
655:   MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
656:   return(0);
657: }

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

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

675: /*
676:    J = [a*I-df_dx, -df_dy
677:         dg_dx, dg_dy]
678: */

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

688:   user->t = t;

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

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

724:   DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet);
725:   DMCompositeScatter(user->dmpgrid,U,Xgen,Xnet);

727:   VecGetArray(Xgen,&xgen);

729:   VecGetArrayRead(U,&u);
730:   VecGetArray(R,&r);
731:   r[0] = 0.;

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

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

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

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

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

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

789:   user.neqs_gen   = 9*ngen; /* # eqs. for generator subsystem */
790:   user.neqs_net   = 2*nbus; /* # eqs. for network subsystem   */
791:   user.neqs_pgrid = user.neqs_gen + user.neqs_net;

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

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

824:   }
825:   PetscOptionsEnd();

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

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

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

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

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

878:   /* SOLVE THE APPLICATION */
879:   TaoSolve(tao); 

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

897: /* ------------------------------------------------------------------ */
898: /*
899:    FormFunction - Evaluates the function and corresponding gradient.

901:    Input Parameters:
902:    tao - the Tao context
903:    X   - the input vector
904:    ptr - optional user-defined context, as set by TaoSetObjectiveAndGradientRoutine()

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

925:   VecGetArrayRead(P,(const PetscScalar**)&x_ptr);
926:   PG[0] = x_ptr[0];
927:   PG[1] = x_ptr[1];
928:   PG[2] = x_ptr[2];
929:   VecRestoreArrayRead(P,(const PetscScalar**)&x_ptr);

931:   ctx->stepnum = 0;

933:   VecZeroEntries(ctx->vec_q);

935:   /* Read initial voltage vector and Ybus */
936:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"X.bin",FILE_MODE_READ,&Xview);
937:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"Ybus.bin",FILE_MODE_READ,&Ybusview);

939:   VecCreate(PETSC_COMM_WORLD,&ctx->V0);
940:   VecSetSizes(ctx->V0,PETSC_DECIDE,ctx->neqs_net);
941:   VecLoad(ctx->V0,Xview);

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

949:   PetscViewerDestroy(&Xview);
950:   PetscViewerDestroy(&Ybusview);

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

954:   MatCreate(PETSC_COMM_WORLD,&J);
955:   MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,ctx->neqs_pgrid,ctx->neqs_pgrid);
956:   MatSetFromOptions(J);
957:   PreallocateJacobian(J,ctx);

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

969:   TSMonitorSet(ts,MonitorUpdateQ,ctx,NULL);
970:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
971:      Set initial conditions
972:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
973:   SetInitialGuess(X,ctx);

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

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

991:   ctx->stepnum++;

993:   TSSetTimeStep(ts,0.01);
994:   TSSetMaxTime(ts,ctx->tmax);
995:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);
996:   TSSetFromOptions(ts);

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

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

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

1022:   MatAssemblyBegin(ctx->Ybus,MAT_FINAL_ASSEMBLY);
1023:   MatAssemblyEnd(ctx->Ybus,MAT_FINAL_ASSEMBLY);

1025:   ctx->alg_flg = PETSC_TRUE;
1026:   /* Solve the algebraic equations */
1027:   SNESSolve(snes_alg,NULL,X);

1029:   ctx->stepnum++;

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

1037:   /* Remove the fault */
1038:   row_loc = 2*ctx->faultbus; col_loc = 2*ctx->faultbus+1;
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;
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:   MatZeroEntries(J);

1050:   ctx->alg_flg = PETSC_TRUE;

1052:   /* Solve the algebraic equations */
1053:   SNESSolve(snes_alg,NULL,X);

1055:   ctx->stepnum++;

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

1063:   VecGetArray(ctx->vec_q,&x_ptr);
1064:   *f   = x_ptr[0];
1065:   VecRestoreArray(ctx->vec_q,&x_ptr);

1067:   MatDestroy(&ctx->Ybus);
1068:   VecDestroy(&ctx->V0);
1069:   SNESDestroy(&snes_alg);
1070:   VecDestroy(&F_alg);
1071:   MatDestroy(&J);
1072:   VecDestroy(&X);
1073:   TSDestroy(&ts);

1075:   return 0;
1076: }

1078: /*TEST

1080:   build:
1081:       requires: double !complex !define(USE_64BIT_INDICES)

1083:    test:
1084:       args: -viewer_binary_skip_info -tao_monitor -tao_gttol .2
1085:       localrunfiles: petscoptions X.bin Ybus.bin

1087: TEST*/