xref: /petsc/src/ts/tutorials/power_grid/ex3.c (revision 6a98f8dc3f2c9149905a87dc2e9d0fedaf64e09a)
1 
2 static char help[] = "Basic equation for generator stability analysis.\n";
3 
4 /*F
5 
6 \begin{eqnarray}
7                  \frac{d \theta}{dt} = \omega_b (\omega - \omega_s)
8                  \frac{2 H}{\omega_s}\frac{d \omega}{dt} & = & P_m - P_max \sin(\theta) -D(\omega - \omega_s)\\
9 \end{eqnarray}
10 
11 
12 
13   Ensemble of initial conditions
14    ./ex3 -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    ./ex3           -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    ./ex3 -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 /*T
35 
36 T*/
37 
38 #include <petscts.h>
39 #include "ex3.h"
40 
41 int main(int argc,char **argv)
42 {
43   TS             ts;            /* ODE integrator */
44   Vec            U;             /* solution will be stored here */
45   Mat            A;             /* Jacobian matrix */
46   PetscErrorCode ierr;
47   PetscMPIInt    size;
48   PetscInt       n = 2;
49   AppCtx         ctx;
50   PetscScalar    *u;
51   PetscReal      du[2] = {0.0,0.0};
52   PetscBool      ensemble = PETSC_FALSE,flg1,flg2;
53   PetscInt       direction[2];
54   PetscBool      terminate[2];
55 
56   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
57      Initialize program
58      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
59   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
60   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
61   if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs");
62 
63   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
64     Create necessary matrix and vectors
65     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
66   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
67   ierr = MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
68   ierr = MatSetType(A,MATDENSE);CHKERRQ(ierr);
69   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
70   ierr = MatSetUp(A);CHKERRQ(ierr);
71 
72   ierr = MatCreateVecs(A,&U,NULL);CHKERRQ(ierr);
73 
74   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
75     Set runtime options
76     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
77   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Swing equation options","");CHKERRQ(ierr);
78   {
79     ctx.omega_b = 1.0;
80     ctx.omega_s = 2.0*PETSC_PI*60.0;
81     ctx.H       = 5.0;
82     ierr        = PetscOptionsScalar("-Inertia","","",ctx.H,&ctx.H,NULL);CHKERRQ(ierr);
83     ctx.D       = 5.0;
84     ierr        = PetscOptionsScalar("-D","","",ctx.D,&ctx.D,NULL);CHKERRQ(ierr);
85     ctx.E       = 1.1378;
86     ctx.V       = 1.0;
87     ctx.X       = 0.545;
88     ctx.Pmax    = ctx.E*ctx.V/ctx.X;
89     ctx.Pmax_ini = ctx.Pmax;
90     ierr        = PetscOptionsScalar("-Pmax","","",ctx.Pmax,&ctx.Pmax,NULL);CHKERRQ(ierr);
91     ctx.Pm      = 0.9;
92     ierr        = PetscOptionsScalar("-Pm","","",ctx.Pm,&ctx.Pm,NULL);CHKERRQ(ierr);
93     ctx.tf      = 1.0;
94     ctx.tcl     = 1.05;
95     ierr        = PetscOptionsReal("-tf","Time to start fault","",ctx.tf,&ctx.tf,NULL);CHKERRQ(ierr);
96     ierr        = PetscOptionsReal("-tcl","Time to end fault","",ctx.tcl,&ctx.tcl,NULL);CHKERRQ(ierr);
97     ierr        = PetscOptionsBool("-ensemble","Run ensemble of different initial conditions","",ensemble,&ensemble,NULL);CHKERRQ(ierr);
98     if (ensemble) {
99       ctx.tf      = -1;
100       ctx.tcl     = -1;
101     }
102 
103     ierr = VecGetArray(U,&u);CHKERRQ(ierr);
104     u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax);
105     u[1] = 1.0;
106     ierr = PetscOptionsRealArray("-u","Initial solution","",u,&n,&flg1);CHKERRQ(ierr);
107     n    = 2;
108     ierr = PetscOptionsRealArray("-du","Perturbation in initial solution","",du,&n,&flg2);CHKERRQ(ierr);
109     u[0] += du[0];
110     u[1] += du[1];
111     ierr = VecRestoreArray(U,&u);CHKERRQ(ierr);
112     if (flg1 || flg2) {
113       ctx.tf      = -1;
114       ctx.tcl     = -1;
115     }
116   }
117   ierr = PetscOptionsEnd();CHKERRQ(ierr);
118 
119   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
120      Create timestepping solver context
121      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
122   ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
123   ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr);
124   ierr = TSSetType(ts,TSTHETA);CHKERRQ(ierr);
125   ierr = TSSetEquationType(ts,TS_EQ_IMPLICIT);CHKERRQ(ierr);
126   ierr = TSARKIMEXSetFullyImplicit(ts,PETSC_TRUE);CHKERRQ(ierr);
127   ierr = TSSetIFunction(ts,NULL,(TSIFunction) IFunction,&ctx);CHKERRQ(ierr);
128   ierr = TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,&ctx);CHKERRQ(ierr);
129 
130   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
131      Set initial conditions
132    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
133   ierr = TSSetSolution(ts,U);CHKERRQ(ierr);
134 
135   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
136      Set solver options
137    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
138   ierr = TSSetMaxTime(ts,35.0);CHKERRQ(ierr);
139   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr);
140   ierr = TSSetTimeStep(ts,.1);CHKERRQ(ierr);
141   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
142 
143   direction[0] = direction[1] = 1;
144   terminate[0] = terminate[1] = PETSC_FALSE;
145 
146   ierr = TSSetEventHandler(ts,2,direction,terminate,EventFunction,PostEventFunction,(void*)&ctx);CHKERRQ(ierr);
147 
148   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
149      Solve nonlinear system
150      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
151   if (ensemble) {
152     for (du[1] = -2.5; du[1] <= .01; du[1] += .1) {
153       ierr = VecGetArray(U,&u);CHKERRQ(ierr);
154       u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax);
155       u[1] = ctx.omega_s;
156       u[0] += du[0];
157       u[1] += du[1];
158       ierr = VecRestoreArray(U,&u);CHKERRQ(ierr);
159       ierr = TSSetTimeStep(ts,.01);CHKERRQ(ierr);
160       ierr = TSSolve(ts,U);CHKERRQ(ierr);
161     }
162   } else {
163     ierr = TSSolve(ts,U);CHKERRQ(ierr);
164   }
165   ierr = VecView(U,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
166   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
167      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
168    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
169   ierr = MatDestroy(&A);CHKERRQ(ierr);
170   ierr = VecDestroy(&U);CHKERRQ(ierr);
171   ierr = TSDestroy(&ts);CHKERRQ(ierr);
172   ierr = PetscFinalize();
173   return ierr;
174 }
175 
176 
177 /*TEST
178 
179    build:
180      requires: !complex !single
181 
182    test:
183       args: -nox
184 
185 TEST*/
186