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