Actual source code: ex2.c


  2: static char help[] = "Basic equation for generator stability analysis.\n";

  4: /*F

  6: \begin{eqnarray}
  7:                  \frac{2 H}{\omega_s}\frac{d \omega}{dt} & = & P_m - \frac{EV}{X} \sin(\theta) -D(\omega - \omega_s)\\
  8:                  \frac{d \theta}{dt} = \omega - \omega_s
  9: \end{eqnarray}

 11:   Ensemble of initial conditions
 12:    ./ex2 -ensemble -ts_monitor_draw_solution_phase -1,-3,3,3      -ts_adapt_dt_max .01  -ts_monitor -ts_type rosw -pc_type lu -ksp_type preonly

 14:   Fault at .1 seconds
 15:    ./ex2           -ts_monitor_draw_solution_phase .42,.95,.6,1.05 -ts_adapt_dt_max .01  -ts_monitor -ts_type rosw -pc_type lu -ksp_type preonly

 17:   Initial conditions same as when fault is ended
 18:    ./ex2 -u 0.496792,1.00932 -ts_monitor_draw_solution_phase .42,.95,.6,1.05  -ts_adapt_dt_max .01  -ts_monitor -ts_type rosw -pc_type lu -ksp_type preonly

 20: F*/

 22: /*
 23:    Include "petscts.h" so that we can use TS solvers.  Note that this
 24:    file automatically includes:
 25:      petscsys.h       - base PETSc routines   petscvec.h - vectors
 26:      petscmat.h - matrices
 27:      petscis.h     - index sets            petscksp.h - Krylov subspace methods
 28:      petscviewer.h - viewers               petscpc.h  - preconditioners
 29:      petscksp.h   - linear solvers
 30: */

 32: #include <petscts.h>

 34: typedef struct {
 35:   PetscScalar H, D, omega_s, Pmax, Pm, E, V, X;
 36:   PetscReal   tf, tcl;
 37: } AppCtx;

 39: /*
 40:      Defines the ODE passed to the ODE solver
 41: */
 42: static PetscErrorCode IFunction(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, AppCtx *ctx)
 43: {
 44:   PetscScalar       *f, Pmax;
 45:   const PetscScalar *u, *udot;

 47:   /*  The next three lines allow us to access the entries of the vectors directly */
 48:   VecGetArrayRead(U, &u);
 49:   VecGetArrayRead(Udot, &udot);
 50:   VecGetArray(F, &f);
 51:   if ((t > ctx->tf) && (t < ctx->tcl)) Pmax = 0.0; /* A short-circuit on the generator terminal that drives the electrical power output (Pmax*sin(delta)) to 0 */
 52:   else if (t >= ctx->tcl) Pmax = ctx->E / 0.745;
 53:   else Pmax = ctx->Pmax;
 54:   f[0] = udot[0] - ctx->omega_s * (u[1] - 1.0);
 55:   f[1] = 2.0 * ctx->H * udot[1] + Pmax * PetscSinScalar(u[0]) + ctx->D * (u[1] - 1.0) - ctx->Pm;

 57:   VecRestoreArrayRead(U, &u);
 58:   VecRestoreArrayRead(Udot, &udot);
 59:   VecRestoreArray(F, &f);
 60:   return 0;
 61: }

 63: /*
 64:      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
 65: */
 66: static PetscErrorCode IJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat B, AppCtx *ctx)
 67: {
 68:   PetscInt           rowcol[] = {0, 1};
 69:   PetscScalar        J[2][2], Pmax;
 70:   const PetscScalar *u, *udot;

 72:   VecGetArrayRead(U, &u);
 73:   VecGetArrayRead(Udot, &udot);
 74:   if ((t > ctx->tf) && (t < ctx->tcl)) Pmax = 0.0; /* A short-circuit on the generator terminal that drives the electrical power output (Pmax*sin(delta)) to 0 */
 75:   else if (t >= ctx->tcl) Pmax = ctx->E / 0.745;
 76:   else Pmax = ctx->Pmax;

 78:   J[0][0] = a;
 79:   J[0][1] = -ctx->omega_s;
 80:   J[1][1] = 2.0 * ctx->H * a + ctx->D;
 81:   J[1][0] = Pmax * PetscCosScalar(u[0]);

 83:   MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES);
 84:   VecRestoreArrayRead(U, &u);
 85:   VecRestoreArrayRead(Udot, &udot);

 87:   MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);
 88:   MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);
 89:   if (A != B) {
 90:     MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY);
 91:     MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY);
 92:   }
 93:   return 0;
 94: }

 96: PetscErrorCode PostStep(TS ts)
 97: {
 98:   Vec       X;
 99:   PetscReal t;

101:   TSGetTime(ts, &t);
102:   if (t >= .2) {
103:     TSGetSolution(ts, &X);
104:     VecView(X, PETSC_VIEWER_STDOUT_WORLD);
105:     exit(0);
106:     /* results in initial conditions after fault of -u 0.496792,1.00932 */
107:   }
108:   return 0;
109: }

111: int main(int argc, char **argv)
112: {
113:   TS           ts; /* ODE integrator */
114:   Vec          U;  /* solution will be stored here */
115:   Mat          A;  /* Jacobian matrix */
116:   PetscMPIInt  size;
117:   PetscInt     n = 2;
118:   AppCtx       ctx;
119:   PetscScalar *u;
120:   PetscReal    du[2]    = {0.0, 0.0};
121:   PetscBool    ensemble = PETSC_FALSE, flg1, flg2;

123:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
124:      Initialize program
125:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
127:   PetscInitialize(&argc, &argv, (char *)0, help);
128:   MPI_Comm_size(PETSC_COMM_WORLD, &size);

131:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
132:     Create necessary matrix and vectors
133:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
134:   MatCreate(PETSC_COMM_WORLD, &A);
135:   MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE);
136:   MatSetType(A, MATDENSE);
137:   MatSetFromOptions(A);
138:   MatSetUp(A);

140:   MatCreateVecs(A, &U, NULL);

142:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
143:     Set runtime options
144:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
145:   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Swing equation options", "");
146:   {
147:     ctx.omega_s = 2.0 * PETSC_PI * 60.0;
148:     ctx.H       = 5.0;
149:     PetscOptionsScalar("-Inertia", "", "", ctx.H, &ctx.H, NULL);
150:     ctx.D = 5.0;
151:     PetscOptionsScalar("-D", "", "", ctx.D, &ctx.D, NULL);
152:     ctx.E    = 1.1378;
153:     ctx.V    = 1.0;
154:     ctx.X    = 0.545;
155:     ctx.Pmax = ctx.E * ctx.V / ctx.X;
156:     PetscOptionsScalar("-Pmax", "", "", ctx.Pmax, &ctx.Pmax, NULL);
157:     ctx.Pm = 0.9;
158:     PetscOptionsScalar("-Pm", "", "", ctx.Pm, &ctx.Pm, NULL);
159:     ctx.tf  = 1.0;
160:     ctx.tcl = 1.05;
161:     PetscOptionsReal("-tf", "Time to start fault", "", ctx.tf, &ctx.tf, NULL);
162:     PetscOptionsReal("-tcl", "Time to end fault", "", ctx.tcl, &ctx.tcl, NULL);
163:     PetscOptionsBool("-ensemble", "Run ensemble of different initial conditions", "", ensemble, &ensemble, NULL);
164:     if (ensemble) {
165:       ctx.tf  = -1;
166:       ctx.tcl = -1;
167:     }

169:     VecGetArray(U, &u);
170:     u[0] = PetscAsinScalar(ctx.Pm / ctx.Pmax);
171:     u[1] = 1.0;
172:     PetscOptionsRealArray("-u", "Initial solution", "", u, &n, &flg1);
173:     n = 2;
174:     PetscOptionsRealArray("-du", "Perturbation in initial solution", "", du, &n, &flg2);
175:     u[0] += du[0];
176:     u[1] += du[1];
177:     VecRestoreArray(U, &u);
178:     if (flg1 || flg2) {
179:       ctx.tf  = -1;
180:       ctx.tcl = -1;
181:     }
182:   }
183:   PetscOptionsEnd();

185:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
186:      Create timestepping solver context
187:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
188:   TSCreate(PETSC_COMM_WORLD, &ts);
189:   TSSetProblemType(ts, TS_NONLINEAR);
190:   TSSetType(ts, TSROSW);
191:   TSSetIFunction(ts, NULL, (TSIFunction)IFunction, &ctx);
192:   TSSetIJacobian(ts, A, A, (TSIJacobian)IJacobian, &ctx);

194:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
195:      Set initial conditions
196:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
197:   TSSetSolution(ts, U);

199:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
200:      Set solver options
201:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
202:   TSSetMaxTime(ts, 35.0);
203:   TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP);
204:   TSSetTimeStep(ts, .01);
205:   TSSetFromOptions(ts);
206:   /* TSSetPostStep(ts,PostStep);  */

208:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
209:      Solve nonlinear system
210:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
211:   if (ensemble) {
212:     for (du[1] = -2.5; du[1] <= .01; du[1] += .1) {
213:       VecGetArray(U, &u);
214:       u[0] = PetscAsinScalar(ctx.Pm / ctx.Pmax);
215:       u[1] = ctx.omega_s;
216:       u[0] += du[0];
217:       u[1] += du[1];
218:       VecRestoreArray(U, &u);
219:       TSSetTimeStep(ts, .01);
220:       TSSolve(ts, U);
221:     }
222:   } else {
223:     TSSolve(ts, U);
224:   }
225:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
226:      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
227:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
228:   MatDestroy(&A);
229:   VecDestroy(&U);
230:   TSDestroy(&ts);
231:   PetscFinalize();
232:   return 0;
233: }

235: /*TEST

237:    build:
238:       requires: !complex

240:    test:
241:       args: -nox -ts_dt 10

243: TEST*/