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