Actual source code: ex2.c


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


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

Ensemble of initial conditions
./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

Fault at .1 seconds
./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

Initial conditions same as when fault is ended
./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

 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;                       J[0][1] = -ctx->omega_s;
 79:   J[1][1] = 2.0*ctx->H*a + ctx->D;   J[1][0] = Pmax*PetscCosScalar(u[0]);

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

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

 94: PetscErrorCode PostStep(TS ts)
 95: {
 96:   Vec            X;
 97:   PetscReal      t;

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

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

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

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

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

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

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

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

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

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

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

233: /*TEST

235:    build:
236:       requires: !complex

238:    test:
239:       args: -nox -ts_dt 10

241: TEST*/