xref: /petsc/src/ts/tutorials/power_grid/ex3.c (revision 8ebe3e4e9e00d86ece2e9fcd0cc84910b0ad437c)
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   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
57   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
58   if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs");
59 
60   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
61     Create necessary matrix and vectors
62     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
63   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
64   ierr = MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
65   ierr = MatSetType(A,MATDENSE);CHKERRQ(ierr);
66   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
67   ierr = MatSetUp(A);CHKERRQ(ierr);
68 
69   ierr = MatCreateVecs(A,&U,NULL);CHKERRQ(ierr);
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     ierr        = PetscOptionsScalar("-Inertia","","",ctx.H,&ctx.H,NULL);CHKERRQ(ierr);
80     ctx.D       = 5.0;
81     ierr        = PetscOptionsScalar("-D","","",ctx.D,&ctx.D,NULL);CHKERRQ(ierr);
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     ierr        = PetscOptionsScalar("-Pmax","","",ctx.Pmax,&ctx.Pmax,NULL);CHKERRQ(ierr);
88     ctx.Pm      = 0.9;
89     ierr        = PetscOptionsScalar("-Pm","","",ctx.Pm,&ctx.Pm,NULL);CHKERRQ(ierr);
90     ctx.tf      = 1.0;
91     ctx.tcl     = 1.05;
92     ierr        = PetscOptionsReal("-tf","Time to start fault","",ctx.tf,&ctx.tf,NULL);CHKERRQ(ierr);
93     ierr        = PetscOptionsReal("-tcl","Time to end fault","",ctx.tcl,&ctx.tcl,NULL);CHKERRQ(ierr);
94     ierr        = PetscOptionsBool("-ensemble","Run ensemble of different initial conditions","",ensemble,&ensemble,NULL);CHKERRQ(ierr);
95     if (ensemble) {
96       ctx.tf      = -1;
97       ctx.tcl     = -1;
98     }
99 
100     ierr = VecGetArray(U,&u);CHKERRQ(ierr);
101     u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax);
102     u[1] = 1.0;
103     ierr = PetscOptionsRealArray("-u","Initial solution","",u,&n,&flg1);CHKERRQ(ierr);
104     n    = 2;
105     ierr = PetscOptionsRealArray("-du","Perturbation in initial solution","",du,&n,&flg2);CHKERRQ(ierr);
106     u[0] += du[0];
107     u[1] += du[1];
108     ierr = VecRestoreArray(U,&u);CHKERRQ(ierr);
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   ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
120   ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr);
121   ierr = TSSetType(ts,TSTHETA);CHKERRQ(ierr);
122   ierr = TSSetEquationType(ts,TS_EQ_IMPLICIT);CHKERRQ(ierr);
123   ierr = TSARKIMEXSetFullyImplicit(ts,PETSC_TRUE);CHKERRQ(ierr);
124   ierr = TSSetIFunction(ts,NULL,(TSIFunction) IFunction,&ctx);CHKERRQ(ierr);
125   ierr = TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,&ctx);CHKERRQ(ierr);
126 
127   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
128      Set initial conditions
129    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
130   ierr = TSSetSolution(ts,U);CHKERRQ(ierr);
131 
132   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
133      Set solver options
134    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
135   ierr = TSSetMaxTime(ts,35.0);CHKERRQ(ierr);
136   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr);
137   ierr = TSSetTimeStep(ts,.1);CHKERRQ(ierr);
138   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
139 
140   direction[0] = direction[1] = 1;
141   terminate[0] = terminate[1] = PETSC_FALSE;
142 
143   ierr = TSSetEventHandler(ts,2,direction,terminate,EventFunction,PostEventFunction,(void*)&ctx);CHKERRQ(ierr);
144 
145   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
146      Solve nonlinear system
147      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
148   if (ensemble) {
149     for (du[1] = -2.5; du[1] <= .01; du[1] += .1) {
150       ierr = VecGetArray(U,&u);CHKERRQ(ierr);
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       ierr = VecRestoreArray(U,&u);CHKERRQ(ierr);
156       ierr = TSSetTimeStep(ts,.01);CHKERRQ(ierr);
157       ierr = TSSolve(ts,U);CHKERRQ(ierr);
158     }
159   } else {
160     ierr = TSSolve(ts,U);CHKERRQ(ierr);
161   }
162   ierr = VecView(U,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
163   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
164      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
165    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
166   ierr = MatDestroy(&A);CHKERRQ(ierr);
167   ierr = VecDestroy(&U);CHKERRQ(ierr);
168   ierr = TSDestroy(&ts);CHKERRQ(ierr);
169   ierr = PetscFinalize();
170   return ierr;
171 }
172 
173 /*TEST
174 
175    build:
176      requires: !complex !single
177 
178    test:
179       args: -nox
180 
181 TEST*/
182