Actual source code: ex9.c
2: static char help[] = "Basic equation for generator stability analysis.\n";
\begin{eqnarray}
\frac{d \theta}{dt} = \omega_b (\omega - \omega_s)
\frac{2 H}{\omega_s}\frac{d \omega}{dt} & = & P_m - P_max \sin(\theta) -D(\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_b,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 RHSFunction(TS ts,PetscReal t,Vec U,Vec F,AppCtx *ctx)
43: {
44: const PetscScalar *u;
45: PetscScalar *f,Pmax;
47: /* The next three lines allow us to access the entries of the vectors directly */
48: VecGetArrayRead(U,&u);
49: VecGetArray(F,&f);
50: 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 */
51: else Pmax = ctx->Pmax;
53: f[0] = ctx->omega_b*(u[1] - ctx->omega_s);
54: f[1] = (-Pmax*PetscSinScalar(u[0]) - ctx->D*(u[1] - ctx->omega_s) + ctx->Pm)*ctx->omega_s/(2.0*ctx->H);
56: VecRestoreArrayRead(U,&u);
57: VecRestoreArray(F,&f);
58: return 0;
59: }
61: /*
62: Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
63: */
64: static PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec U,Mat A,Mat B,AppCtx *ctx)
65: {
66: PetscInt rowcol[] = {0,1};
67: PetscScalar J[2][2],Pmax;
68: const PetscScalar *u;
70: VecGetArrayRead(U,&u);
71: 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 */
72: else Pmax = ctx->Pmax;
74: J[0][0] = 0; J[0][1] = ctx->omega_b;
75: J[1][1] = -ctx->D*ctx->omega_s/(2.0*ctx->H); J[1][0] = -Pmax*PetscCosScalar(u[0])*ctx->omega_s/(2.0*ctx->H);
77: MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);
78: VecRestoreArrayRead(U,&u);
80: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
81: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
82: if (A != B) {
83: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
84: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
85: }
86: return 0;
87: }
89: int main(int argc,char **argv)
90: {
91: TS ts; /* ODE integrator */
92: Vec U; /* solution will be stored here */
93: Mat A; /* Jacobian matrix */
95: PetscMPIInt size;
96: PetscInt n = 2;
97: AppCtx ctx;
98: PetscScalar *u;
99: PetscReal du[2] = {0.0,0.0};
100: PetscBool ensemble = PETSC_FALSE,flg1,flg2;
102: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
103: Initialize program
104: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
105: PetscInitialize(&argc,&argv,(char*)0,help);
106: MPI_Comm_size(PETSC_COMM_WORLD,&size);
109: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
110: Create necessary matrix and vectors
111: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
112: MatCreate(PETSC_COMM_WORLD,&A);
113: MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);
114: MatSetType(A,MATDENSE);
115: MatSetFromOptions(A);
116: MatSetUp(A);
118: MatCreateVecs(A,&U,NULL);
120: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
121: Set runtime options
122: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
123: PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Swing equation options","");
124: {
125: ctx.omega_b = 1.0;
126: ctx.omega_s = 2.0*PETSC_PI*60.0;
127: ctx.H = 5.0;
128: PetscOptionsScalar("-Inertia","","",ctx.H,&ctx.H,NULL);
129: ctx.D = 5.0;
130: PetscOptionsScalar("-D","","",ctx.D,&ctx.D,NULL);
131: ctx.E = 1.1378;
132: ctx.V = 1.0;
133: ctx.X = 0.545;
134: ctx.Pmax = ctx.E*ctx.V/ctx.X;
135: PetscOptionsScalar("-Pmax","","",ctx.Pmax,&ctx.Pmax,NULL);
136: ctx.Pm = 0.9;
137: PetscOptionsScalar("-Pm","","",ctx.Pm,&ctx.Pm,NULL);
138: ctx.tf = 1.0;
139: ctx.tcl = 1.05;
140: PetscOptionsReal("-tf","Time to start fault","",ctx.tf,&ctx.tf,NULL);
141: PetscOptionsReal("-tcl","Time to end fault","",ctx.tcl,&ctx.tcl,NULL);
142: PetscOptionsBool("-ensemble","Run ensemble of different initial conditions","",ensemble,&ensemble,NULL);
143: if (ensemble) {
144: ctx.tf = -1;
145: ctx.tcl = -1;
146: }
148: VecGetArray(U,&u);
149: u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax);
150: u[1] = 1.0;
151: PetscOptionsRealArray("-u","Initial solution","",u,&n,&flg1);
152: n = 2;
153: PetscOptionsRealArray("-du","Perturbation in initial solution","",du,&n,&flg2);
154: u[0] += du[0];
155: u[1] += du[1];
156: VecRestoreArray(U,&u);
157: if (flg1 || flg2) {
158: ctx.tf = -1;
159: ctx.tcl = -1;
160: }
161: }
162: PetscOptionsEnd();
164: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
165: Create timestepping solver context
166: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
167: TSCreate(PETSC_COMM_WORLD,&ts);
168: TSSetProblemType(ts,TS_NONLINEAR);
169: TSSetType(ts,TSTHETA);
170: TSSetRHSFunction(ts,NULL,(TSRHSFunction)RHSFunction,&ctx);
171: TSSetRHSJacobian(ts,A,A,(TSRHSJacobian)RHSJacobian,&ctx);
173: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
174: Set initial conditions
175: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
176: TSSetSolution(ts,U);
178: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
179: Set solver options
180: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
181: TSSetMaxTime(ts,35.0);
182: TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);
183: TSSetTimeStep(ts,.01);
184: TSSetFromOptions(ts);
186: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
187: Solve nonlinear system
188: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
189: if (ensemble) {
190: for (du[1] = -2.5; du[1] <= .01; du[1] += .1) {
191: VecGetArray(U,&u);
192: u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax);
193: u[1] = ctx.omega_s;
194: u[0] += du[0];
195: u[1] += du[1];
196: VecRestoreArray(U,&u);
197: TSSetTimeStep(ts,.01);
198: TSSolve(ts,U);
199: }
200: } else {
201: TSSolve(ts,U);
202: }
203: VecView(U,PETSC_VIEWER_STDOUT_WORLD);
204: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
205: Free work space. All PETSc objects should be destroyed when they are no longer needed.
206: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
207: MatDestroy(&A);
208: VecDestroy(&U);
209: TSDestroy(&ts);
210: PetscFinalize();
211: return 0;
212: }
214: /*TEST
216: build:
217: requires: !complex
219: test:
221: TEST*/