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: PetscErrorCode ierr;
45: const PetscScalar *u;
46: PetscScalar *f,Pmax;
49: /* The next three lines allow us to access the entries of the vectors directly */
50: VecGetArrayRead(U,&u);
51: VecGetArray(F,&f);
52: 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 */
53: else Pmax = ctx->Pmax;
55: f[0] = ctx->omega_b*(u[1] - ctx->omega_s);
56: f[1] = (-Pmax*PetscSinScalar(u[0]) - ctx->D*(u[1] - ctx->omega_s) + ctx->Pm)*ctx->omega_s/(2.0*ctx->H);
58: VecRestoreArrayRead(U,&u);
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 RHSJacobian(TS ts,PetscReal t,Vec U,Mat A,Mat B,AppCtx *ctx)
67: {
68: PetscErrorCode ierr;
69: PetscInt rowcol[] = {0,1};
70: PetscScalar J[2][2],Pmax;
71: const PetscScalar *u;
74: VecGetArrayRead(U,&u);
75: 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 */
76: else Pmax = ctx->Pmax;
78: J[0][0] = 0; J[0][1] = ctx->omega_b;
79: 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);
81: MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);
82: VecRestoreArrayRead(U,&u);
84: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
85: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
86: if (A != B) {
87: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
88: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
89: }
90: return(0);
91: }
93: int main(int argc,char **argv)
94: {
95: TS ts; /* ODE integrator */
96: Vec U; /* solution will be stored here */
97: Mat A; /* Jacobian matrix */
99: PetscMPIInt size;
100: PetscInt n = 2;
101: AppCtx ctx;
102: PetscScalar *u;
103: PetscReal du[2] = {0.0,0.0};
104: PetscBool ensemble = PETSC_FALSE,flg1,flg2;
106: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
107: Initialize program
108: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
109: PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
110: MPI_Comm_size(PETSC_COMM_WORLD,&size);
111: if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs");
113: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
114: Create necessary matrix and vectors
115: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
116: MatCreate(PETSC_COMM_WORLD,&A);
117: MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);
118: MatSetType(A,MATDENSE);
119: MatSetFromOptions(A);
120: MatSetUp(A);
122: MatCreateVecs(A,&U,NULL);
124: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
125: Set runtime options
126: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
127: PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Swing equation options","");
128: {
129: ctx.omega_b = 1.0;
130: ctx.omega_s = 2.0*PETSC_PI*60.0;
131: ctx.H = 5.0;
132: PetscOptionsScalar("-Inertia","","",ctx.H,&ctx.H,NULL);
133: ctx.D = 5.0;
134: PetscOptionsScalar("-D","","",ctx.D,&ctx.D,NULL);
135: ctx.E = 1.1378;
136: ctx.V = 1.0;
137: ctx.X = 0.545;
138: ctx.Pmax = ctx.E*ctx.V/ctx.X;
139: PetscOptionsScalar("-Pmax","","",ctx.Pmax,&ctx.Pmax,NULL);
140: ctx.Pm = 0.9;
141: PetscOptionsScalar("-Pm","","",ctx.Pm,&ctx.Pm,NULL);
142: ctx.tf = 1.0;
143: ctx.tcl = 1.05;
144: PetscOptionsReal("-tf","Time to start fault","",ctx.tf,&ctx.tf,NULL);
145: PetscOptionsReal("-tcl","Time to end fault","",ctx.tcl,&ctx.tcl,NULL);
146: PetscOptionsBool("-ensemble","Run ensemble of different initial conditions","",ensemble,&ensemble,NULL);
147: if (ensemble) {
148: ctx.tf = -1;
149: ctx.tcl = -1;
150: }
152: VecGetArray(U,&u);
153: u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax);
154: u[1] = 1.0;
155: PetscOptionsRealArray("-u","Initial solution","",u,&n,&flg1);
156: n = 2;
157: PetscOptionsRealArray("-du","Perturbation in initial solution","",du,&n,&flg2);
158: u[0] += du[0];
159: u[1] += du[1];
160: VecRestoreArray(U,&u);
161: if (flg1 || flg2) {
162: ctx.tf = -1;
163: ctx.tcl = -1;
164: }
165: }
166: PetscOptionsEnd();
168: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
169: Create timestepping solver context
170: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
171: TSCreate(PETSC_COMM_WORLD,&ts);
172: TSSetProblemType(ts,TS_NONLINEAR);
173: TSSetType(ts,TSTHETA);
174: TSSetRHSFunction(ts,NULL,(TSRHSFunction)RHSFunction,&ctx);
175: TSSetRHSJacobian(ts,A,A,(TSRHSJacobian)RHSJacobian,&ctx);
177: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
178: Set initial conditions
179: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
180: TSSetSolution(ts,U);
182: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
183: Set solver options
184: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
185: TSSetMaxTime(ts,35.0);
186: TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);
187: TSSetTimeStep(ts,.01);
188: TSSetFromOptions(ts);
190: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
191: Solve nonlinear system
192: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
193: if (ensemble) {
194: for (du[1] = -2.5; du[1] <= .01; du[1] += .1) {
195: VecGetArray(U,&u);
196: u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax);
197: u[1] = ctx.omega_s;
198: u[0] += du[0];
199: u[1] += du[1];
200: VecRestoreArray(U,&u);
201: TSSetTimeStep(ts,.01);
202: TSSolve(ts,U);
203: }
204: } else {
205: TSSolve(ts,U);
206: }
207: VecView(U,PETSC_VIEWER_STDOUT_WORLD);
208: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
209: Free work space. All PETSc objects should be destroyed when they are no longer needed.
210: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
211: MatDestroy(&A);
212: VecDestroy(&U);
213: TSDestroy(&ts);
214: PetscFinalize();
215: return ierr;
216: }
218: /*TEST
220: build:
221: requires: !complex
223: test:
225: TEST*/