Actual source code: ex3.c
petsc-3.9.4 2018-09-11
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: /*T
36: T*/
38: #include <petscts.h>
40: typedef struct {
41: PetscScalar H,D,omega_b,omega_s,Pmax,Pmax_ini,Pm,E,V,X;
42: PetscReal tf,tcl;
43: } AppCtx;
45: /* Event check */
46: PetscErrorCode EventFunction(TS ts,PetscReal t,Vec X,PetscScalar *fvalue,void *ctx)
47: {
48: AppCtx *user=(AppCtx*)ctx;
51: /* Event for fault-on time */
52: fvalue[0] = t - user->tf;
53: /* Event for fault-off time */
54: fvalue[1] = t - user->tcl;
56: return (0);
57: }
59: PetscErrorCode PostEventFunction(TS ts,PetscInt nevents,PetscInt event_list[],PetscReal t,Vec X,PetscBool forwardsolve,void* ctx)
60: {
61: AppCtx *user=(AppCtx*)ctx;
62:
65: if (event_list[0] == 0) user->Pmax = 0.0; /* Apply disturbance - this is done by setting Pmax = 0 */
66: else if (event_list[0] == 1) user->Pmax = user->Pmax_ini; /* Remove the fault - this is done by setting Pmax = Pmax_ini */
67: return (0);
68: }
70: /*
71: Defines the ODE passed to the ODE solver
72: */
73: static PetscErrorCode IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,AppCtx *ctx)
74: {
75: PetscErrorCode ierr;
76: const PetscScalar *u,*udot;
77: PetscScalar *f,Pmax;
80: /* The next three lines allow us to access the entries of the vectors directly */
81: VecGetArrayRead (U,&u);
82: VecGetArrayRead (Udot,&udot);
83: VecGetArray (F,&f);
84: Pmax = ctx->Pmax;
86: f[0] = udot[0] - ctx->omega_b*(u[1] - ctx->omega_s);
87: f[1] = 2.0*ctx->H/ctx->omega_s*udot[1] + Pmax*PetscSinScalar(u[0]) + ctx->D*(u[1] - ctx->omega_s)- ctx->Pm;
89: VecRestoreArrayRead (U,&u);
90: VecRestoreArrayRead (Udot,&udot);
91: VecRestoreArray (F,&f);
92: return (0);
93: }
95: /*
96: Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian () for the meaning of a and the Jacobian.
97: */
98: static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,AppCtx *ctx)
99: {
100: PetscErrorCode ierr;
101: PetscInt rowcol[] = {0,1};
102: PetscScalar J[2][2],Pmax;
103: const PetscScalar *u,*udot;
106: VecGetArrayRead (U,&u);
107: VecGetArrayRead (Udot,&udot);
108: Pmax = ctx->Pmax;
110: J[0][0] = a; J[0][1] = -ctx->omega_b;
111: J[1][1] = 2.0*ctx->H/ctx->omega_s*a + ctx->D; J[1][0] = Pmax*PetscCosScalar(u[0]);
113: MatSetValues (B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES );
114: VecRestoreArrayRead (U,&u);
115: VecRestoreArrayRead (Udot,&udot);
117: MatAssemblyBegin (A,MAT_FINAL_ASSEMBLY );
118: MatAssemblyEnd (A,MAT_FINAL_ASSEMBLY );
119: if (A != B) {
120: MatAssemblyBegin (B,MAT_FINAL_ASSEMBLY );
121: MatAssemblyEnd (B,MAT_FINAL_ASSEMBLY );
122: }
123: return (0);
124: }
126: int main(int argc,char **argv)
127: {
128: TS ts; /* ODE integrator */
129: Vec U; /* solution will be stored here */
130: Mat A; /* Jacobian matrix */
132: PetscMPIInt size;
133: PetscInt n = 2;
134: AppCtx ctx;
135: PetscScalar *u;
136: PetscReal du[2] = {0.0,0.0};
137: PetscBool ensemble = PETSC_FALSE ,flg1,flg2;
138: PetscInt direction[2];
139: PetscBool terminate[2];
141: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
142: Initialize program
143: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
144: PetscInitialize (&argc,&argv,(char*)0,help);if (ierr) return ierr;
145: MPI_Comm_size (PETSC_COMM_WORLD ,&size);
146: if (size > 1) SETERRQ (PETSC_COMM_WORLD ,PETSC_ERR_SUP,"Only for sequential runs" );
148: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
149: Create necessary matrix and vectors
150: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
151: MatCreate (PETSC_COMM_WORLD ,&A);
152: MatSetSizes (A,n,n,PETSC_DETERMINE ,PETSC_DETERMINE );
153: MatSetType (A,MATDENSE );
154: MatSetFromOptions (A);
155: MatSetUp (A);
157: MatCreateVecs (A,&U,NULL);
159: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
160: Set runtime options
161: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
162: PetscOptionsBegin (PETSC_COMM_WORLD ,NULL,"Swing equation options" ,"" );
163: {
164: ctx.omega_b = 1.0;
165: ctx.omega_s = 2.0*PETSC_PI*60.0;
166: ctx.H = 5.0;
167: PetscOptionsScalar ("-Inertia" ,"" ,"" ,ctx.H,&ctx.H,NULL);
168: ctx.D = 5.0;
169: PetscOptionsScalar ("-D" ,"" ,"" ,ctx.D,&ctx.D,NULL);
170: ctx.E = 1.1378;
171: ctx.V = 1.0;
172: ctx.X = 0.545;
173: ctx.Pmax = ctx.E*ctx.V/ctx.X;
174: ctx.Pmax_ini = ctx.Pmax;
175: PetscOptionsScalar ("-Pmax" ,"" ,"" ,ctx.Pmax,&ctx.Pmax,NULL);
176: ctx.Pm = 0.9;
177: PetscOptionsScalar ("-Pm" ,"" ,"" ,ctx.Pm,&ctx.Pm,NULL);
178: ctx.tf = 1.0;
179: ctx.tcl = 1.05;
180: PetscOptionsReal ("-tf" ,"Time to start fault" ,"" ,ctx.tf,&ctx.tf,NULL);
181: PetscOptionsReal ("-tcl" ,"Time to end fault" ,"" ,ctx.tcl,&ctx.tcl,NULL);
182: PetscOptionsBool ("-ensemble" ,"Run ensemble of different initial conditions" ,"" ,ensemble,&ensemble,NULL);
183: if (ensemble) {
184: ctx.tf = -1;
185: ctx.tcl = -1;
186: }
188: VecGetArray (U,&u);
189: u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax);
190: u[1] = 1.0;
191: PetscOptionsRealArray ("-u" ,"Initial solution" ,"" ,u,&n,&flg1);
192: n = 2;
193: PetscOptionsRealArray ("-du" ,"Perturbation in initial solution" ,"" ,du,&n,&flg2);
194: u[0] += du[0];
195: u[1] += du[1];
196: VecRestoreArray (U,&u);
197: if (flg1 || flg2) {
198: ctx.tf = -1;
199: ctx.tcl = -1;
200: }
201: }
202: PetscOptionsEnd ();
204: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
205: Create timestepping solver context
206: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
207: TSCreate (PETSC_COMM_WORLD ,&ts);
208: TSSetProblemType (ts,TS_NONLINEAR );
209: TSSetType (ts,TSTHETA );
210: TSSetEquationType (ts,TS_EQ_IMPLICIT );
211: TSARKIMEXSetFullyImplicit (ts,PETSC_TRUE );
212: TSSetIFunction (ts,NULL,(TSIFunction) IFunction,&ctx);
213: TSSetIJacobian (ts,A,A,(TSIJacobian)IJacobian,&ctx);
215: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
216: Set initial conditions
217: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
218: TSSetSolution (ts,U);
220: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
221: Set solver options
222: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
223: TSSetMaxTime (ts,35.0);
224: TSSetExactFinalTime (ts,TS_EXACTFINALTIME_MATCHSTEP );
225: TSSetTimeStep (ts,.1);
226: TSSetFromOptions (ts);
228: direction[0] = direction[1] = 1;
229: terminate[0] = terminate[1] = PETSC_FALSE ;
231: TSSetEventHandler (ts,2,direction,terminate,EventFunction,PostEventFunction,(void*)&ctx);
233: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
234: Solve nonlinear system
235: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
236: if (ensemble) {
237: for (du[1] = -2.5; du[1] <= .01; du[1] += .1) {
238: VecGetArray (U,&u);
239: u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax);
240: u[1] = ctx.omega_s;
241: u[0] += du[0];
242: u[1] += du[1];
243: VecRestoreArray (U,&u);
244: TSSetTimeStep (ts,.01);
245: TSSolve (ts,U);
246: }
247: } else {
248: TSSolve (ts,U);
249: }
250: VecView (U,PETSC_VIEWER_STDOUT_WORLD );
251: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
252: Free work space. All PETSc objects should be destroyed when they are no longer needed.
253: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
254: MatDestroy (&A);
255: VecDestroy (&U);
256: TSDestroy (&ts);
257: PetscFinalize ();
258: return ierr;
259: }
262: /*TEST
264: build:
265: requires: !complex !single
267: test:
268: args: -nox
270: TEST*/