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