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