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