Actual source code: ex2.c
petsc-3.12.1 2019-10-22
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
25: /*
26: Include "petscts.h" so that we can use TS solvers. Note that this
27: file automatically includes:
28: petscsys.h - base PETSc routines petscvec.h - vectors
29: petscmat.h - matrices
30: petscis.h - index sets petscksp.h - Krylov subspace methods
31: petscviewer.h - viewers petscpc.h - preconditioners
32: petscksp.h - linear solvers
33: */
35: #include <petscts.h>
37: typedef struct {
38: PetscScalar H,D,omega_s,Pmax,Pm,E,V,X;
39: PetscReal tf,tcl;
40: } AppCtx;
42: /*
43: Defines the ODE passed to the ODE solver
44: */
45: static PetscErrorCode IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,AppCtx *ctx)
46: {
47: PetscErrorCode ierr;
48: PetscScalar *f,Pmax;
49: const PetscScalar *u,*udot;
52: /* The next three lines allow us to access the entries of the vectors directly */
53: VecGetArrayRead (U,&u);
54: VecGetArrayRead (Udot,&udot);
55: VecGetArray (F,&f);
56: 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 */
57: else if (t >= ctx->tcl) Pmax = ctx->E/0.745;
58: else Pmax = ctx->Pmax;
59: f[0] = udot[0] - ctx->omega_s*(u[1] - 1.0);
60: f[1] = 2.0*ctx->H*udot[1] + Pmax*PetscSinScalar(u[0]) + ctx->D*(u[1] - 1.0)- ctx->Pm;
62: VecRestoreArrayRead (U,&u);
63: VecRestoreArrayRead (Udot,&udot);
64: VecRestoreArray (F,&f);
65: return (0);
66: }
68: /*
69: Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian () for the meaning of a and the Jacobian.
70: */
71: static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,AppCtx *ctx)
72: {
73: PetscErrorCode ierr;
74: PetscInt rowcol[] = {0,1};
75: PetscScalar J[2][2],Pmax;
76: const PetscScalar *u,*udot;
79: VecGetArrayRead (U,&u);
80: VecGetArrayRead (Udot,&udot);
81: 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 */
82: else if (t >= ctx->tcl) Pmax = ctx->E/0.745;
83: else Pmax = ctx->Pmax;
85: J[0][0] = a; J[0][1] = -ctx->omega_s;
86: J[1][1] = 2.0*ctx->H*a + ctx->D; J[1][0] = Pmax*PetscCosScalar(u[0]);
88: MatSetValues (B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES );
89: VecRestoreArrayRead (U,&u);
90: VecRestoreArrayRead (Udot,&udot);
92: MatAssemblyBegin (A,MAT_FINAL_ASSEMBLY );
93: MatAssemblyEnd (A,MAT_FINAL_ASSEMBLY );
94: if (A != B) {
95: MatAssemblyBegin (B,MAT_FINAL_ASSEMBLY );
96: MatAssemblyEnd (B,MAT_FINAL_ASSEMBLY );
97: }
98: return (0);
99: }
102: PetscErrorCode PostStep(TS ts)
103: {
105: Vec X;
106: PetscReal t;
109: TSGetTime (ts,&t);
110: if (t >= .2) {
111: TSGetSolution (ts,&X);
112: VecView (X,PETSC_VIEWER_STDOUT_WORLD );
113: exit(0);
114: /* results in initial conditions after fault of -u 0.496792,1.00932 */
115: }
116: return (0);
117: }
120: int main(int argc,char **argv)
121: {
122: TS ts; /* ODE integrator */
123: Vec U; /* solution will be stored here */
124: Mat A; /* Jacobian matrix */
126: PetscMPIInt size;
127: PetscInt n = 2;
128: AppCtx ctx;
129: PetscScalar *u;
130: PetscReal du[2] = {0.0,0.0};
131: PetscBool ensemble = PETSC_FALSE ,flg1,flg2;
133: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
134: Initialize program
135: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
136: PetscInitialize (&argc,&argv,(char*)0,help);if (ierr) return ierr;
137: MPI_Comm_size (PETSC_COMM_WORLD ,&size);
138: if (size > 1) SETERRQ (PETSC_COMM_WORLD ,PETSC_ERR_SUP,"Only for sequential runs" );
140: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
141: Create necessary matrix and vectors
142: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
143: MatCreate (PETSC_COMM_WORLD ,&A);
144: MatSetSizes (A,n,n,PETSC_DETERMINE ,PETSC_DETERMINE );
145: MatSetType (A,MATDENSE );
146: MatSetFromOptions (A);
147: MatSetUp (A);
149: MatCreateVecs (A,&U,NULL);
151: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
152: Set runtime options
153: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
154: PetscOptionsBegin (PETSC_COMM_WORLD ,NULL,"Swing equation options" ,"" );
155: {
156: ctx.omega_s = 2.0*PETSC_PI*60.0;
157: ctx.H = 5.0;
158: PetscOptionsScalar ("-Inertia" ,"" ,"" ,ctx.H,&ctx.H,NULL);
159: ctx.D = 5.0;
160: PetscOptionsScalar ("-D" ,"" ,"" ,ctx.D,&ctx.D,NULL);
161: ctx.E = 1.1378;
162: ctx.V = 1.0;
163: ctx.X = 0.545;
164: ctx.Pmax = ctx.E*ctx.V/ctx.X;
165: PetscOptionsScalar ("-Pmax" ,"" ,"" ,ctx.Pmax,&ctx.Pmax,NULL);
166: ctx.Pm = 0.9;
167: PetscOptionsScalar ("-Pm" ,"" ,"" ,ctx.Pm,&ctx.Pm,NULL);
168: ctx.tf = 1.0;
169: ctx.tcl = 1.05;
170: PetscOptionsReal ("-tf" ,"Time to start fault" ,"" ,ctx.tf,&ctx.tf,NULL);
171: PetscOptionsReal ("-tcl" ,"Time to end fault" ,"" ,ctx.tcl,&ctx.tcl,NULL);
172: PetscOptionsBool ("-ensemble" ,"Run ensemble of different initial conditions" ,"" ,ensemble,&ensemble,NULL);
173: if (ensemble) {
174: ctx.tf = -1;
175: ctx.tcl = -1;
176: }
178: VecGetArray (U,&u);
179: u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax);
180: u[1] = 1.0;
181: PetscOptionsRealArray ("-u" ,"Initial solution" ,"" ,u,&n,&flg1);
182: n = 2;
183: PetscOptionsRealArray ("-du" ,"Perturbation in initial solution" ,"" ,du,&n,&flg2);
184: u[0] += du[0];
185: u[1] += du[1];
186: VecRestoreArray (U,&u);
187: if (flg1 || flg2) {
188: ctx.tf = -1;
189: ctx.tcl = -1;
190: }
191: }
192: PetscOptionsEnd ();
194: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
195: Create timestepping solver context
196: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
197: TSCreate (PETSC_COMM_WORLD ,&ts);
198: TSSetProblemType (ts,TS_NONLINEAR );
199: TSSetType (ts,TSROSW );
200: TSSetIFunction (ts,NULL,(TSIFunction) IFunction,&ctx);
201: TSSetIJacobian (ts,A,A,(TSIJacobian)IJacobian,&ctx);
203: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
204: Set initial conditions
205: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
206: TSSetSolution (ts,U);
208: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
209: Set solver options
210: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
211: TSSetMaxTime (ts,35.0);
212: TSSetExactFinalTime (ts,TS_EXACTFINALTIME_MATCHSTEP );
213: TSSetTimeStep (ts,.01);
214: TSSetFromOptions (ts);
215: /* TSSetPostStep (ts,PostStep); */
218: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
219: Solve nonlinear system
220: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
221: if (ensemble) {
222: for (du[1] = -2.5; du[1] <= .01; du[1] += .1) {
223: VecGetArray (U,&u);
224: u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax);
225: u[1] = ctx.omega_s;
226: u[0] += du[0];
227: u[1] += du[1];
228: VecRestoreArray (U,&u);
229: TSSetTimeStep (ts,.01);
230: TSSolve (ts,U);
231: }
232: } else {
233: TSSolve (ts,U);
234: }
235: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
236: Free work space. All PETSc objects should be destroyed when they are no longer needed.
237: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
238: MatDestroy (&A);
239: VecDestroy (&U);
240: TSDestroy (&ts);
241: PetscFinalize ();
242: return ierr;
243: }
246: /*TEST
248: build:
249: requires: !complex
251: test:
252: args: -nox -ts_dt 10
254: TEST*/