xref: /petsc/src/ts/tutorials/power_grid/ex1.c (revision dbbe0bcd3f3a8fbab5a45420dc06f8387e5764c6)
1 
2 static char help[] = "Basic equation for generator stability analysis.\n";
3 
4 /*F
5 
6 \begin{eqnarray}
7                  \frac{2 H}{\omega_s}\frac{d \omega}{dt} & = & P_m - \frac{EV}{X} \sin(\theta) \\
8                  \frac{d \theta}{dt} = \omega - \omega_s
9 \end{eqnarray}
10 
11 F*/
12 
13 /*
14    Include "petscts.h" so that we can use TS solvers.  Note that this
15    file automatically includes:
16      petscsys.h       - base PETSc routines   petscvec.h - vectors
17      petscmat.h - matrices
18      petscis.h     - index sets            petscksp.h - Krylov subspace methods
19      petscviewer.h - viewers               petscpc.h  - preconditioners
20      petscksp.h   - linear solvers
21 */
22 
23 #include <petscts.h>
24 
25 typedef struct {
26   PetscScalar H,omega_s,E,V,X;
27   PetscRandom rand;
28 } AppCtx;
29 
30 /*
31      Defines the ODE passed to the ODE solver
32 */
33 static PetscErrorCode IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,AppCtx *ctx)
34 {
35   PetscScalar        *f,r;
36   const PetscScalar  *u,*udot;
37   static PetscScalar R = .4;
38 
39   PetscFunctionBegin;
40   PetscCall(PetscRandomGetValue(ctx->rand,&r));
41   if (r > .9) R = .5;
42   if (r < .1) R = .4;
43   R = .4;
44   /*  The next three lines allow us to access the entries of the vectors directly */
45   PetscCall(VecGetArrayRead(U,&u));
46   PetscCall(VecGetArrayRead(Udot,&udot));
47   PetscCall(VecGetArray(F,&f));
48   f[0] = 2.0*ctx->H*udot[0]/ctx->omega_s + ctx->E*ctx->V*PetscSinScalar(u[1])/ctx->X - R;
49   f[1] = udot[1] - u[0] + ctx->omega_s;
50 
51   PetscCall(VecRestoreArrayRead(U,&u));
52   PetscCall(VecRestoreArrayRead(Udot,&udot));
53   PetscCall(VecRestoreArray(F,&f));
54   PetscFunctionReturn(0);
55 }
56 
57 /*
58      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
59 */
60 static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,AppCtx *ctx)
61 {
62   PetscInt          rowcol[] = {0,1};
63   PetscScalar       J[2][2];
64   const PetscScalar *u,*udot;
65 
66   PetscFunctionBegin;
67   PetscCall(VecGetArrayRead(U,&u));
68   PetscCall(VecGetArrayRead(Udot,&udot));
69   J[0][0] = 2.0*ctx->H*a/ctx->omega_s;   J[0][1] = -ctx->E*ctx->V*PetscCosScalar(u[1])/ctx->X;
70   J[1][0] = -1.0;                        J[1][1] = a;
71   PetscCall(MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES));
72   PetscCall(VecRestoreArrayRead(U,&u));
73   PetscCall(VecRestoreArrayRead(Udot,&udot));
74 
75   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
76   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
77   if (A != B) {
78     PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
79     PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
80   }
81   PetscFunctionReturn(0);
82 }
83 
84 int main(int argc,char **argv)
85 {
86   TS             ts;            /* ODE integrator */
87   Vec            U;             /* solution will be stored here */
88   Mat            A;             /* Jacobian matrix */
89   PetscMPIInt    size;
90   PetscInt       n = 2;
91   AppCtx         ctx;
92   PetscScalar    *u;
93 
94   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
95      Initialize program
96      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
97   PetscFunctionBeginUser;
98   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
99   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
100   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Only for sequential runs");
101 
102   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
103     Create necessary matrix and vectors
104     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
105   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
106   PetscCall(MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE));
107   PetscCall(MatSetFromOptions(A));
108   PetscCall(MatSetUp(A));
109 
110   PetscCall(MatCreateVecs(A,&U,NULL));
111 
112   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
113     Set runtime options
114     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
115   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Reaction options","");
116   {
117     ctx.omega_s = 1.0;
118     PetscCall(PetscOptionsScalar("-omega_s","","",ctx.omega_s,&ctx.omega_s,NULL));
119     ctx.H       = 1.0;
120     PetscCall(PetscOptionsScalar("-H","","",ctx.H,&ctx.H,NULL));
121     ctx.E       = 1.0;
122     PetscCall(PetscOptionsScalar("-E","","",ctx.E,&ctx.E,NULL));
123     ctx.V       = 1.0;
124     PetscCall(PetscOptionsScalar("-V","","",ctx.V,&ctx.V,NULL));
125     ctx.X       = 1.0;
126     PetscCall(PetscOptionsScalar("-X","","",ctx.X,&ctx.X,NULL));
127 
128     PetscCall(VecGetArray(U,&u));
129     u[0] = 1;
130     u[1] = .7;
131     PetscCall(VecRestoreArray(U,&u));
132     PetscCall(PetscOptionsGetVec(NULL,NULL,"-initial",U,NULL));
133   }
134   PetscOptionsEnd();
135 
136   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&ctx.rand));
137   PetscCall(PetscRandomSetFromOptions(ctx.rand));
138 
139   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
140      Create timestepping solver context
141      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
142   PetscCall(TSCreate(PETSC_COMM_WORLD,&ts));
143   PetscCall(TSSetProblemType(ts,TS_NONLINEAR));
144   PetscCall(TSSetType(ts,TSROSW));
145   PetscCall(TSSetIFunction(ts,NULL,(TSIFunction) IFunction,&ctx));
146   PetscCall(TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,&ctx));
147 
148   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
149      Set initial conditions
150    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
151   PetscCall(TSSetSolution(ts,U));
152 
153   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
154      Set solver options
155    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
156   PetscCall(TSSetMaxTime(ts,2000.0));
157   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP));
158   PetscCall(TSSetTimeStep(ts,.001));
159   PetscCall(TSSetFromOptions(ts));
160 
161   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
162      Solve nonlinear system
163      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
164   PetscCall(TSSolve(ts,U));
165 
166   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
167      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
168    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
169   PetscCall(MatDestroy(&A));
170   PetscCall(VecDestroy(&U));
171   PetscCall(TSDestroy(&ts));
172   PetscCall(PetscRandomDestroy(&ctx.rand));
173   PetscCall(PetscFinalize());
174   return 0;
175 }
176 
177 /*TEST
178 
179    build:
180      requires: !complex !single
181 
182    test:
183       args: -ksp_gmres_cgs_refinement_type refine_always -snes_type newtonls -ts_max_steps 10
184 
185    test:
186       suffix: 2
187       args: -ts_max_steps 10
188 
189 TEST*/
190