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