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