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