xref: /petsc/src/ts/tutorials/power_grid/ex9.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    ./ex2 -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    ./ex2           -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    ./ex2 -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 
35*c4762a1bSJed Brown #include <petscts.h>
36*c4762a1bSJed Brown 
37*c4762a1bSJed Brown typedef struct {
38*c4762a1bSJed Brown   PetscScalar H,D,omega_b,omega_s,Pmax,Pm,E,V,X;
39*c4762a1bSJed Brown   PetscReal   tf,tcl;
40*c4762a1bSJed Brown } AppCtx;
41*c4762a1bSJed Brown 
42*c4762a1bSJed Brown /*
43*c4762a1bSJed Brown      Defines the ODE passed to the ODE solver
44*c4762a1bSJed Brown */
45*c4762a1bSJed Brown static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec U,Vec F,AppCtx *ctx)
46*c4762a1bSJed Brown {
47*c4762a1bSJed Brown   PetscErrorCode    ierr;
48*c4762a1bSJed Brown   const PetscScalar *u;
49*c4762a1bSJed Brown   PetscScalar       *f,Pmax;
50*c4762a1bSJed Brown 
51*c4762a1bSJed Brown   PetscFunctionBegin;
52*c4762a1bSJed Brown   /*  The next three lines allow us to access the entries of the vectors directly */
53*c4762a1bSJed Brown   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
54*c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
55*c4762a1bSJed Brown   if ((t > ctx->tf) && (t < ctx->tcl)) Pmax = 0.0; /* A short-circuit on the generator terminal that drives the electrical power output (Pmax*sin(delta)) to 0 */
56*c4762a1bSJed Brown   else Pmax = ctx->Pmax;
57*c4762a1bSJed Brown 
58*c4762a1bSJed Brown   f[0] = ctx->omega_b*(u[1] - ctx->omega_s);
59*c4762a1bSJed Brown   f[1] = (-Pmax*PetscSinScalar(u[0]) - ctx->D*(u[1] - ctx->omega_s) + ctx->Pm)*ctx->omega_s/(2.0*ctx->H);
60*c4762a1bSJed Brown 
61*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
62*c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
63*c4762a1bSJed Brown   PetscFunctionReturn(0);
64*c4762a1bSJed Brown }
65*c4762a1bSJed Brown 
66*c4762a1bSJed Brown /*
67*c4762a1bSJed Brown      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
68*c4762a1bSJed Brown */
69*c4762a1bSJed Brown static PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec U,Mat A,Mat B,AppCtx *ctx)
70*c4762a1bSJed Brown {
71*c4762a1bSJed Brown   PetscErrorCode    ierr;
72*c4762a1bSJed Brown   PetscInt          rowcol[] = {0,1};
73*c4762a1bSJed Brown   PetscScalar       J[2][2],Pmax;
74*c4762a1bSJed Brown   const PetscScalar *u;
75*c4762a1bSJed Brown 
76*c4762a1bSJed Brown   PetscFunctionBegin;
77*c4762a1bSJed Brown   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
78*c4762a1bSJed Brown   if ((t > ctx->tf) && (t < ctx->tcl)) Pmax = 0.0; /* A short-circuit on the generator terminal that drives the electrical power output (Pmax*sin(delta)) to 0 */
79*c4762a1bSJed Brown   else Pmax = ctx->Pmax;
80*c4762a1bSJed Brown 
81*c4762a1bSJed Brown   J[0][0] = 0;                                    J[0][1] = ctx->omega_b;
82*c4762a1bSJed Brown   J[1][1] = -ctx->D*ctx->omega_s/(2.0*ctx->H);    J[1][0] = -Pmax*PetscCosScalar(u[0])*ctx->omega_s/(2.0*ctx->H);
83*c4762a1bSJed Brown 
84*c4762a1bSJed Brown   ierr    = MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);CHKERRQ(ierr);
85*c4762a1bSJed Brown   ierr    = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
86*c4762a1bSJed Brown 
87*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
88*c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
89*c4762a1bSJed Brown   if (A != B) {
90*c4762a1bSJed Brown     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
91*c4762a1bSJed Brown     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
92*c4762a1bSJed Brown   }
93*c4762a1bSJed Brown   PetscFunctionReturn(0);
94*c4762a1bSJed Brown }
95*c4762a1bSJed Brown 
96*c4762a1bSJed Brown int main(int argc,char **argv)
97*c4762a1bSJed Brown {
98*c4762a1bSJed Brown   TS             ts;            /* ODE integrator */
99*c4762a1bSJed Brown   Vec            U;             /* solution will be stored here */
100*c4762a1bSJed Brown   Mat            A;             /* Jacobian matrix */
101*c4762a1bSJed Brown   PetscErrorCode ierr;
102*c4762a1bSJed Brown   PetscMPIInt    size;
103*c4762a1bSJed Brown   PetscInt       n = 2;
104*c4762a1bSJed Brown   AppCtx         ctx;
105*c4762a1bSJed Brown   PetscScalar    *u;
106*c4762a1bSJed Brown   PetscReal      du[2] = {0.0,0.0};
107*c4762a1bSJed Brown   PetscBool      ensemble = PETSC_FALSE,flg1,flg2;
108*c4762a1bSJed Brown 
109*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
110*c4762a1bSJed Brown      Initialize program
111*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
112*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
113*c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
114*c4762a1bSJed Brown   if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs");
115*c4762a1bSJed Brown 
116*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
117*c4762a1bSJed Brown     Create necessary matrix and vectors
118*c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
119*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
120*c4762a1bSJed Brown   ierr = MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
121*c4762a1bSJed Brown   ierr = MatSetType(A,MATDENSE);CHKERRQ(ierr);
122*c4762a1bSJed Brown   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
123*c4762a1bSJed Brown   ierr = MatSetUp(A);CHKERRQ(ierr);
124*c4762a1bSJed Brown 
125*c4762a1bSJed Brown   ierr = MatCreateVecs(A,&U,NULL);CHKERRQ(ierr);
126*c4762a1bSJed Brown 
127*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
128*c4762a1bSJed Brown     Set runtime options
129*c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
130*c4762a1bSJed Brown   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Swing equation options","");CHKERRQ(ierr);
131*c4762a1bSJed Brown   {
132*c4762a1bSJed Brown     ctx.omega_b = 1.0;
133*c4762a1bSJed Brown     ctx.omega_s = 2.0*PETSC_PI*60.0;
134*c4762a1bSJed Brown     ctx.H       = 5.0;
135*c4762a1bSJed Brown     ierr        = PetscOptionsScalar("-Inertia","","",ctx.H,&ctx.H,NULL);CHKERRQ(ierr);
136*c4762a1bSJed Brown     ctx.D       = 5.0;
137*c4762a1bSJed Brown     ierr        = PetscOptionsScalar("-D","","",ctx.D,&ctx.D,NULL);CHKERRQ(ierr);
138*c4762a1bSJed Brown     ctx.E       = 1.1378;
139*c4762a1bSJed Brown     ctx.V       = 1.0;
140*c4762a1bSJed Brown     ctx.X       = 0.545;
141*c4762a1bSJed Brown     ctx.Pmax    = ctx.E*ctx.V/ctx.X;
142*c4762a1bSJed Brown     ierr        = PetscOptionsScalar("-Pmax","","",ctx.Pmax,&ctx.Pmax,NULL);CHKERRQ(ierr);
143*c4762a1bSJed Brown     ctx.Pm      = 0.9;
144*c4762a1bSJed Brown     ierr        = PetscOptionsScalar("-Pm","","",ctx.Pm,&ctx.Pm,NULL);CHKERRQ(ierr);
145*c4762a1bSJed Brown     ctx.tf      = 1.0;
146*c4762a1bSJed Brown     ctx.tcl     = 1.05;
147*c4762a1bSJed Brown     ierr        = PetscOptionsReal("-tf","Time to start fault","",ctx.tf,&ctx.tf,NULL);CHKERRQ(ierr);
148*c4762a1bSJed Brown     ierr        = PetscOptionsReal("-tcl","Time to end fault","",ctx.tcl,&ctx.tcl,NULL);CHKERRQ(ierr);
149*c4762a1bSJed Brown     ierr        = PetscOptionsBool("-ensemble","Run ensemble of different initial conditions","",ensemble,&ensemble,NULL);CHKERRQ(ierr);
150*c4762a1bSJed Brown     if (ensemble) {
151*c4762a1bSJed Brown       ctx.tf      = -1;
152*c4762a1bSJed Brown       ctx.tcl     = -1;
153*c4762a1bSJed Brown     }
154*c4762a1bSJed Brown 
155*c4762a1bSJed Brown     ierr = VecGetArray(U,&u);CHKERRQ(ierr);
156*c4762a1bSJed Brown     u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax);
157*c4762a1bSJed Brown     u[1] = 1.0;
158*c4762a1bSJed Brown     ierr = PetscOptionsRealArray("-u","Initial solution","",u,&n,&flg1);CHKERRQ(ierr);
159*c4762a1bSJed Brown     n    = 2;
160*c4762a1bSJed Brown     ierr = PetscOptionsRealArray("-du","Perturbation in initial solution","",du,&n,&flg2);CHKERRQ(ierr);
161*c4762a1bSJed Brown     u[0] += du[0];
162*c4762a1bSJed Brown     u[1] += du[1];
163*c4762a1bSJed Brown     ierr = VecRestoreArray(U,&u);CHKERRQ(ierr);
164*c4762a1bSJed Brown     if (flg1 || flg2) {
165*c4762a1bSJed Brown       ctx.tf      = -1;
166*c4762a1bSJed Brown       ctx.tcl     = -1;
167*c4762a1bSJed Brown     }
168*c4762a1bSJed Brown   }
169*c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
170*c4762a1bSJed Brown 
171*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
172*c4762a1bSJed Brown      Create timestepping solver context
173*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
174*c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
175*c4762a1bSJed Brown   ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr);
176*c4762a1bSJed Brown   ierr = TSSetType(ts,TSTHETA);CHKERRQ(ierr);
177*c4762a1bSJed Brown   ierr = TSSetRHSFunction(ts,NULL,(TSRHSFunction)RHSFunction,&ctx);CHKERRQ(ierr);
178*c4762a1bSJed Brown   ierr = TSSetRHSJacobian(ts,A,A,(TSRHSJacobian)RHSJacobian,&ctx);CHKERRQ(ierr);
179*c4762a1bSJed Brown 
180*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
181*c4762a1bSJed Brown      Set initial conditions
182*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
183*c4762a1bSJed Brown   ierr = TSSetSolution(ts,U);CHKERRQ(ierr);
184*c4762a1bSJed Brown 
185*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
186*c4762a1bSJed Brown      Set solver options
187*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
188*c4762a1bSJed Brown   ierr = TSSetMaxTime(ts,35.0);CHKERRQ(ierr);
189*c4762a1bSJed Brown   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr);
190*c4762a1bSJed Brown   ierr = TSSetTimeStep(ts,.01);CHKERRQ(ierr);
191*c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
192*c4762a1bSJed Brown 
193*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
194*c4762a1bSJed Brown      Solve nonlinear system
195*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
196*c4762a1bSJed Brown   if (ensemble) {
197*c4762a1bSJed Brown     for (du[1] = -2.5; du[1] <= .01; du[1] += .1) {
198*c4762a1bSJed Brown       ierr = VecGetArray(U,&u);CHKERRQ(ierr);
199*c4762a1bSJed Brown       u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax);
200*c4762a1bSJed Brown       u[1] = ctx.omega_s;
201*c4762a1bSJed Brown       u[0] += du[0];
202*c4762a1bSJed Brown       u[1] += du[1];
203*c4762a1bSJed Brown       ierr = VecRestoreArray(U,&u);CHKERRQ(ierr);
204*c4762a1bSJed Brown       ierr = TSSetTimeStep(ts,.01);CHKERRQ(ierr);
205*c4762a1bSJed Brown       ierr = TSSolve(ts,U);CHKERRQ(ierr);
206*c4762a1bSJed Brown     }
207*c4762a1bSJed Brown   } else {
208*c4762a1bSJed Brown     ierr = TSSolve(ts,U);CHKERRQ(ierr);
209*c4762a1bSJed Brown   }
210*c4762a1bSJed Brown   ierr = VecView(U,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
211*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
212*c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
213*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
214*c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
215*c4762a1bSJed Brown   ierr = VecDestroy(&U);CHKERRQ(ierr);
216*c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
217*c4762a1bSJed Brown   ierr = PetscFinalize();
218*c4762a1bSJed Brown   return ierr;
219*c4762a1bSJed Brown }
220*c4762a1bSJed Brown 
221*c4762a1bSJed Brown 
222*c4762a1bSJed Brown /*TEST
223*c4762a1bSJed Brown 
224*c4762a1bSJed Brown    build:
225*c4762a1bSJed Brown      requires: !complex
226*c4762a1bSJed Brown 
227*c4762a1bSJed Brown    test:
228*c4762a1bSJed Brown 
229*c4762a1bSJed Brown TEST*/
230