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:   PetscErrorCode    ierr;
 45:   PetscScalar       *f,Pmax;
 46:   const PetscScalar *u,*udot;

 49:   /*  The next three lines allow us to access the entries of the vectors directly */
 50:   VecGetArrayRead(U,&u);
 51:   VecGetArrayRead(Udot,&udot);
 52:   VecGetArray(F,&f);
 53:   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 */
 54:   else if (t >= ctx->tcl) Pmax = ctx->E/0.745;
 55:   else Pmax = ctx->Pmax;
 56:   f[0] = udot[0] - ctx->omega_s*(u[1] - 1.0);
 57:   f[1] = 2.0*ctx->H*udot[1] +  Pmax*PetscSinScalar(u[0]) + ctx->D*(u[1] - 1.0)- ctx->Pm;

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

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

 76:   VecGetArrayRead(U,&u);
 77:   VecGetArrayRead(Udot,&udot);
 78:   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 */
 79:   else if (t >= ctx->tcl) Pmax = ctx->E/0.745;
 80:   else Pmax = ctx->Pmax;

 82:   J[0][0] = a;                       J[0][1] = -ctx->omega_s;
 83:   J[1][1] = 2.0*ctx->H*a + ctx->D;   J[1][0] = Pmax*PetscCosScalar(u[0]);

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

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

 98: PetscErrorCode PostStep(TS ts)
 99: {
101:   Vec            X;
102:   PetscReal      t;

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

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

128:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
129:      Initialize program
130:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
131:   PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
132:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
133:   if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs");

135:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
136:     Create necessary matrix and vectors
137:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
138:   MatCreate(PETSC_COMM_WORLD,&A);
139:   MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);
140:   MatSetType(A,MATDENSE);
141:   MatSetFromOptions(A);
142:   MatSetUp(A);

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

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

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

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

198:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
199:      Set initial conditions
200:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
201:   TSSetSolution(ts,U);

203:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
204:      Set solver options
205:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
206:   TSSetMaxTime(ts,35.0);
207:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);
208:   TSSetTimeStep(ts,.01);
209:   TSSetFromOptions(ts);
210:   /* TSSetPostStep(ts,PostStep);  */

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

239: /*TEST

241:    build:
242:       requires: !complex

244:    test:
245:       args: -nox -ts_dt 10

247: TEST*/