xref: /petsc/src/ts/tutorials/power_grid/ex1.c (revision 6a98f8dc3f2c9149905a87dc2e9d0fedaf64e09a)
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   PetscErrorCode     ierr;
36   PetscScalar        *f,r;
37   const PetscScalar  *u,*udot;
38   static PetscScalar R = .4;
39 
40   PetscFunctionBegin;
41   ierr = PetscRandomGetValue(ctx->rand,&r);CHKERRQ(ierr);
42   if (r > .9) R = .5;
43   if (r < .1) R = .4;
44   R = .4;
45   /*  The next three lines allow us to access the entries of the vectors directly */
46   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
47   ierr = VecGetArrayRead(Udot,&udot);CHKERRQ(ierr);
48   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
49   f[0] = 2.0*ctx->H*udot[0]/ctx->omega_s + ctx->E*ctx->V*PetscSinScalar(u[1])/ctx->X - R;
50   f[1] = udot[1] - u[0] + ctx->omega_s;
51 
52   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
53   ierr = VecRestoreArrayRead(Udot,&udot);CHKERRQ(ierr);
54   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
55   PetscFunctionReturn(0);
56 }
57 
58 /*
59      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
60 */
61 static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,AppCtx *ctx)
62 {
63   PetscErrorCode    ierr;
64   PetscInt          rowcol[] = {0,1};
65   PetscScalar       J[2][2];
66   const PetscScalar *u,*udot;
67 
68   PetscFunctionBegin;
69   ierr    = VecGetArrayRead(U,&u);CHKERRQ(ierr);
70   ierr    = VecGetArrayRead(Udot,&udot);CHKERRQ(ierr);
71   J[0][0] = 2.0*ctx->H*a/ctx->omega_s;   J[0][1] = -ctx->E*ctx->V*PetscCosScalar(u[1])/ctx->X;
72   J[1][0] = -1.0;                        J[1][1] = a;
73   ierr    = MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);CHKERRQ(ierr);
74   ierr    = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
75   ierr    = VecRestoreArrayRead(Udot,&udot);CHKERRQ(ierr);
76 
77   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
78   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
79   if (A != B) {
80     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
81     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
82   }
83   PetscFunctionReturn(0);
84 }
85 
86 int main(int argc,char **argv)
87 {
88   TS             ts;            /* ODE integrator */
89   Vec            U;             /* solution will be stored here */
90   Mat            A;             /* Jacobian matrix */
91   PetscErrorCode ierr;
92   PetscMPIInt    size;
93   PetscInt       n = 2;
94   AppCtx         ctx;
95   PetscScalar    *u;
96 
97   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
98      Initialize program
99      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
100   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
101   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
102   if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs");
103 
104   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
105     Create necessary matrix and vectors
106     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
107   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
108   ierr = MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
109   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
110   ierr = MatSetUp(A);CHKERRQ(ierr);
111 
112   ierr = MatCreateVecs(A,&U,NULL);CHKERRQ(ierr);
113 
114   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
115     Set runtime options
116     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
117   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Reaction options","");CHKERRQ(ierr);
118   {
119     ctx.omega_s = 1.0;
120     ierr        = PetscOptionsScalar("-omega_s","","",ctx.omega_s,&ctx.omega_s,NULL);CHKERRQ(ierr);
121     ctx.H       = 1.0;
122     ierr        = PetscOptionsScalar("-H","","",ctx.H,&ctx.H,NULL);CHKERRQ(ierr);
123     ctx.E       = 1.0;
124     ierr        = PetscOptionsScalar("-E","","",ctx.E,&ctx.E,NULL);CHKERRQ(ierr);
125     ctx.V       = 1.0;
126     ierr        = PetscOptionsScalar("-V","","",ctx.V,&ctx.V,NULL);CHKERRQ(ierr);
127     ctx.X       = 1.0;
128     ierr        = PetscOptionsScalar("-X","","",ctx.X,&ctx.X,NULL);CHKERRQ(ierr);
129 
130     ierr = VecGetArray(U,&u);CHKERRQ(ierr);
131     u[0] = 1;
132     u[1] = .7;
133     ierr = VecRestoreArray(U,&u);CHKERRQ(ierr);
134     ierr = PetscOptionsGetVec(NULL,NULL,"-initial",U,NULL);CHKERRQ(ierr);
135   }
136   ierr = PetscOptionsEnd();CHKERRQ(ierr);
137 
138   ierr = PetscRandomCreate(PETSC_COMM_WORLD,&ctx.rand);CHKERRQ(ierr);
139   ierr = PetscRandomSetFromOptions(ctx.rand);CHKERRQ(ierr);
140 
141   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
142      Create timestepping solver context
143      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
144   ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
145   ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr);
146   ierr = TSSetType(ts,TSROSW);CHKERRQ(ierr);
147   ierr = TSSetIFunction(ts,NULL,(TSIFunction) IFunction,&ctx);CHKERRQ(ierr);
148   ierr = TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,&ctx);CHKERRQ(ierr);
149 
150   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
151      Set initial conditions
152    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
153   ierr = TSSetSolution(ts,U);CHKERRQ(ierr);
154 
155   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
156      Set solver options
157    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
158   ierr = TSSetMaxTime(ts,2000.0);CHKERRQ(ierr);
159   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr);
160   ierr = TSSetTimeStep(ts,.001);CHKERRQ(ierr);
161   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
162 
163   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
164      Solve nonlinear system
165      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
166   ierr = TSSolve(ts,U);CHKERRQ(ierr);
167 
168   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
169      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
170    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
171   ierr = MatDestroy(&A);CHKERRQ(ierr);
172   ierr = VecDestroy(&U);CHKERRQ(ierr);
173   ierr = TSDestroy(&ts);CHKERRQ(ierr);
174   ierr = PetscRandomDestroy(&ctx.rand);CHKERRQ(ierr);
175   ierr = PetscFinalize();
176   return ierr;
177 }
178 
179 
180 /*TEST
181 
182    build:
183      requires: !complex !single
184 
185    test:
186       args: -ksp_gmres_cgs_refinement_type refine_always -snes_type newtonls -ts_max_steps 10
187 
188    test:
189       suffix: 2
190       args: -ts_max_steps 10
191 
192 TEST*/
193