xref: /petsc/src/ts/tutorials/power_grid/ex2.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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{2 H}{\omega_s}\frac{d \omega}{dt} & = & P_m - \frac{EV}{X} \sin(\theta) -D(\omega - \omega_s)\\
8c4762a1bSJed Brown                  \frac{d \theta}{dt} = \omega - \omega_s
9c4762a1bSJed Brown \end{eqnarray}
10c4762a1bSJed Brown 
11c4762a1bSJed Brown   Ensemble of initial conditions
12c4762a1bSJed 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
13c4762a1bSJed Brown 
14c4762a1bSJed Brown   Fault at .1 seconds
15c4762a1bSJed 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
16c4762a1bSJed Brown 
17c4762a1bSJed Brown   Initial conditions same as when fault is ended
18c4762a1bSJed 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
19c4762a1bSJed Brown 
20c4762a1bSJed Brown F*/
21c4762a1bSJed Brown 
22c4762a1bSJed Brown /*
23c4762a1bSJed Brown    Include "petscts.h" so that we can use TS solvers.  Note that this
24c4762a1bSJed Brown    file automatically includes:
25c4762a1bSJed Brown      petscsys.h       - base PETSc routines   petscvec.h - vectors
26c4762a1bSJed Brown      petscmat.h - matrices
27c4762a1bSJed Brown      petscis.h     - index sets            petscksp.h - Krylov subspace methods
28c4762a1bSJed Brown      petscviewer.h - viewers               petscpc.h  - preconditioners
29c4762a1bSJed Brown      petscksp.h   - linear solvers
30c4762a1bSJed Brown */
31c4762a1bSJed Brown 
32c4762a1bSJed Brown #include <petscts.h>
33c4762a1bSJed Brown 
34c4762a1bSJed Brown typedef struct {
35c4762a1bSJed Brown   PetscScalar H,D,omega_s,Pmax,Pm,E,V,X;
36c4762a1bSJed Brown   PetscReal   tf,tcl;
37c4762a1bSJed Brown } AppCtx;
38c4762a1bSJed Brown 
39c4762a1bSJed Brown /*
40c4762a1bSJed Brown      Defines the ODE passed to the ODE solver
41c4762a1bSJed Brown */
42c4762a1bSJed Brown static PetscErrorCode IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,AppCtx *ctx)
43c4762a1bSJed Brown {
44c4762a1bSJed Brown   PetscScalar       *f,Pmax;
45c4762a1bSJed Brown   const PetscScalar *u,*udot;
46c4762a1bSJed Brown 
47c4762a1bSJed Brown   PetscFunctionBegin;
48c4762a1bSJed Brown   /*  The next three lines allow us to access the entries of the vectors directly */
49*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(U,&u));
50*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(Udot,&udot));
51*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(F,&f));
52c4762a1bSJed 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 */
53c4762a1bSJed Brown   else if (t >= ctx->tcl) Pmax = ctx->E/0.745;
54c4762a1bSJed Brown   else Pmax = ctx->Pmax;
55c4762a1bSJed Brown   f[0] = udot[0] - ctx->omega_s*(u[1] - 1.0);
56c4762a1bSJed Brown   f[1] = 2.0*ctx->H*udot[1] +  Pmax*PetscSinScalar(u[0]) + ctx->D*(u[1] - 1.0)- ctx->Pm;
57c4762a1bSJed Brown 
58*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(U,&u));
59*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(Udot,&udot));
60*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(F,&f));
61c4762a1bSJed Brown   PetscFunctionReturn(0);
62c4762a1bSJed Brown }
63c4762a1bSJed Brown 
64c4762a1bSJed Brown /*
65c4762a1bSJed Brown      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
66c4762a1bSJed Brown */
67c4762a1bSJed Brown static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,AppCtx *ctx)
68c4762a1bSJed Brown {
69c4762a1bSJed Brown   PetscInt          rowcol[] = {0,1};
70c4762a1bSJed Brown   PetscScalar       J[2][2],Pmax;
71c4762a1bSJed Brown   const PetscScalar *u,*udot;
72c4762a1bSJed Brown 
73c4762a1bSJed Brown   PetscFunctionBegin;
74*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(U,&u));
75*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(Udot,&udot));
76c4762a1bSJed 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 */
77c4762a1bSJed Brown   else if (t >= ctx->tcl) Pmax = ctx->E/0.745;
78c4762a1bSJed Brown   else Pmax = ctx->Pmax;
79c4762a1bSJed Brown 
80c4762a1bSJed Brown   J[0][0] = a;                       J[0][1] = -ctx->omega_s;
81c4762a1bSJed Brown   J[1][1] = 2.0*ctx->H*a + ctx->D;   J[1][0] = Pmax*PetscCosScalar(u[0]);
82c4762a1bSJed Brown 
83*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES));
84*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(U,&u));
85*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(Udot,&udot));
86c4762a1bSJed Brown 
87*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
88*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
89c4762a1bSJed Brown   if (A != B) {
90*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
91*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
92c4762a1bSJed Brown   }
93c4762a1bSJed Brown   PetscFunctionReturn(0);
94c4762a1bSJed Brown }
95c4762a1bSJed Brown 
96c4762a1bSJed Brown PetscErrorCode PostStep(TS ts)
97c4762a1bSJed Brown {
98c4762a1bSJed Brown   Vec            X;
99c4762a1bSJed Brown   PetscReal      t;
100c4762a1bSJed Brown 
101c4762a1bSJed Brown   PetscFunctionBegin;
102*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetTime(ts,&t));
103c4762a1bSJed Brown   if (t >= .2) {
104*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSGetSolution(ts,&X));
105*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecView(X,PETSC_VIEWER_STDOUT_WORLD));
106c4762a1bSJed Brown     exit(0);
107c4762a1bSJed Brown     /* results in initial conditions after fault of -u 0.496792,1.00932 */
108c4762a1bSJed Brown   }
109c4762a1bSJed Brown   PetscFunctionReturn(0);
110c4762a1bSJed Brown }
111c4762a1bSJed Brown 
112c4762a1bSJed Brown int main(int argc,char **argv)
113c4762a1bSJed Brown {
114c4762a1bSJed Brown   TS             ts;            /* ODE integrator */
115c4762a1bSJed Brown   Vec            U;             /* solution will be stored here */
116c4762a1bSJed Brown   Mat            A;             /* Jacobian matrix */
117c4762a1bSJed Brown   PetscErrorCode ierr;
118c4762a1bSJed Brown   PetscMPIInt    size;
119c4762a1bSJed Brown   PetscInt       n = 2;
120c4762a1bSJed Brown   AppCtx         ctx;
121c4762a1bSJed Brown   PetscScalar    *u;
122c4762a1bSJed Brown   PetscReal      du[2] = {0.0,0.0};
123c4762a1bSJed Brown   PetscBool      ensemble = PETSC_FALSE,flg1,flg2;
124c4762a1bSJed Brown 
125c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
126c4762a1bSJed Brown      Initialize program
127c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
128c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
129*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
1303c633725SBarry Smith   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Only for sequential runs");
131c4762a1bSJed Brown 
132c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
133c4762a1bSJed Brown     Create necessary matrix and vectors
134c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
135*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
136*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE));
137*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(A,MATDENSE));
138*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(A));
139*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(A));
140c4762a1bSJed Brown 
141*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(A,&U,NULL));
142c4762a1bSJed Brown 
143c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
144c4762a1bSJed Brown     Set runtime options
145c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
146c4762a1bSJed Brown   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Swing equation options","");CHKERRQ(ierr);
147c4762a1bSJed Brown   {
148c4762a1bSJed Brown     ctx.omega_s = 2.0*PETSC_PI*60.0;
149c4762a1bSJed Brown     ctx.H       = 5.0;
150*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsScalar("-Inertia","","",ctx.H,&ctx.H,NULL));
151c4762a1bSJed Brown     ctx.D       = 5.0;
152*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsScalar("-D","","",ctx.D,&ctx.D,NULL));
153c4762a1bSJed Brown     ctx.E       = 1.1378;
154c4762a1bSJed Brown     ctx.V       = 1.0;
155c4762a1bSJed Brown     ctx.X       = 0.545;
156c4762a1bSJed Brown     ctx.Pmax    = ctx.E*ctx.V/ctx.X;
157*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsScalar("-Pmax","","",ctx.Pmax,&ctx.Pmax,NULL));
158c4762a1bSJed Brown     ctx.Pm      = 0.9;
159*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsScalar("-Pm","","",ctx.Pm,&ctx.Pm,NULL));
160c4762a1bSJed Brown     ctx.tf      = 1.0;
161c4762a1bSJed Brown     ctx.tcl     = 1.05;
162*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsReal("-tf","Time to start fault","",ctx.tf,&ctx.tf,NULL));
163*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsReal("-tcl","Time to end fault","",ctx.tcl,&ctx.tcl,NULL));
164*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsBool("-ensemble","Run ensemble of different initial conditions","",ensemble,&ensemble,NULL));
165c4762a1bSJed Brown     if (ensemble) {
166c4762a1bSJed Brown       ctx.tf      = -1;
167c4762a1bSJed Brown       ctx.tcl     = -1;
168c4762a1bSJed Brown     }
169c4762a1bSJed Brown 
170*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetArray(U,&u));
171c4762a1bSJed Brown     u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax);
172c4762a1bSJed Brown     u[1] = 1.0;
173*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsRealArray("-u","Initial solution","",u,&n,&flg1));
174c4762a1bSJed Brown     n    = 2;
175*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsRealArray("-du","Perturbation in initial solution","",du,&n,&flg2));
176c4762a1bSJed Brown     u[0] += du[0];
177c4762a1bSJed Brown     u[1] += du[1];
178*5f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreArray(U,&u));
179c4762a1bSJed Brown     if (flg1 || flg2) {
180c4762a1bSJed Brown       ctx.tf      = -1;
181c4762a1bSJed Brown       ctx.tcl     = -1;
182c4762a1bSJed Brown     }
183c4762a1bSJed Brown   }
184c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
185c4762a1bSJed Brown 
186c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
187c4762a1bSJed Brown      Create timestepping solver context
188c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
189*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts));
190*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetProblemType(ts,TS_NONLINEAR));
191*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetType(ts,TSROSW));
192*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetIFunction(ts,NULL,(TSIFunction) IFunction,&ctx));
193*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,&ctx));
194c4762a1bSJed Brown 
195c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
196c4762a1bSJed Brown      Set initial conditions
197c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
198*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetSolution(ts,U));
199c4762a1bSJed Brown 
200c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
201c4762a1bSJed Brown      Set solver options
202c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
203*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetMaxTime(ts,35.0));
204*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP));
205*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTimeStep(ts,.01));
206*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(ts));
207*5f80ce2aSJacob Faibussowitsch   /* CHKERRQ(TSSetPostStep(ts,PostStep));  */
208c4762a1bSJed Brown 
209c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
210c4762a1bSJed Brown      Solve nonlinear system
211c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
212c4762a1bSJed Brown   if (ensemble) {
213c4762a1bSJed Brown     for (du[1] = -2.5; du[1] <= .01; du[1] += .1) {
214*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecGetArray(U,&u));
215c4762a1bSJed Brown       u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax);
216c4762a1bSJed Brown       u[1] = ctx.omega_s;
217c4762a1bSJed Brown       u[0] += du[0];
218c4762a1bSJed Brown       u[1] += du[1];
219*5f80ce2aSJacob Faibussowitsch       CHKERRQ(VecRestoreArray(U,&u));
220*5f80ce2aSJacob Faibussowitsch       CHKERRQ(TSSetTimeStep(ts,.01));
221*5f80ce2aSJacob Faibussowitsch       CHKERRQ(TSSolve(ts,U));
222c4762a1bSJed Brown     }
223c4762a1bSJed Brown   } else {
224*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSolve(ts,U));
225c4762a1bSJed Brown   }
226c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
227c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
228c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
229*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
230*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&U));
231*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDestroy(&ts));
232c4762a1bSJed Brown   ierr = PetscFinalize();
233c4762a1bSJed Brown   return ierr;
234c4762a1bSJed Brown }
235c4762a1bSJed Brown 
236c4762a1bSJed Brown /*TEST
237c4762a1bSJed Brown 
238c4762a1bSJed Brown    build:
239c4762a1bSJed Brown       requires: !complex
240c4762a1bSJed Brown 
241c4762a1bSJed Brown    test:
242c4762a1bSJed Brown       args: -nox -ts_dt 10
243c4762a1bSJed Brown 
244c4762a1bSJed Brown TEST*/
245