xref: /petsc/src/ts/tutorials/power_grid/ex1.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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) \\
8c4762a1bSJed Brown                  \frac{d \theta}{dt} = \omega - \omega_s
9c4762a1bSJed Brown \end{eqnarray}
10c4762a1bSJed Brown 
11c4762a1bSJed Brown F*/
12c4762a1bSJed Brown 
13c4762a1bSJed Brown /*
14c4762a1bSJed Brown    Include "petscts.h" so that we can use TS solvers.  Note that this
15c4762a1bSJed Brown    file automatically includes:
16c4762a1bSJed Brown      petscsys.h       - base PETSc routines   petscvec.h - vectors
17c4762a1bSJed Brown      petscmat.h - matrices
18c4762a1bSJed Brown      petscis.h     - index sets            petscksp.h - Krylov subspace methods
19c4762a1bSJed Brown      petscviewer.h - viewers               petscpc.h  - preconditioners
20c4762a1bSJed Brown      petscksp.h   - linear solvers
21c4762a1bSJed Brown */
22c4762a1bSJed Brown 
23c4762a1bSJed Brown #include <petscts.h>
24c4762a1bSJed Brown 
25c4762a1bSJed Brown typedef struct {
26c4762a1bSJed Brown   PetscScalar H,omega_s,E,V,X;
27c4762a1bSJed Brown   PetscRandom rand;
28c4762a1bSJed Brown } AppCtx;
29c4762a1bSJed Brown 
30c4762a1bSJed Brown /*
31c4762a1bSJed Brown      Defines the ODE passed to the ODE solver
32c4762a1bSJed Brown */
33c4762a1bSJed Brown static PetscErrorCode IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,AppCtx *ctx)
34c4762a1bSJed Brown {
35c4762a1bSJed Brown   PetscScalar        *f,r;
36c4762a1bSJed Brown   const PetscScalar  *u,*udot;
37c4762a1bSJed Brown   static PetscScalar R = .4;
38c4762a1bSJed Brown 
39c4762a1bSJed Brown   PetscFunctionBegin;
405f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomGetValue(ctx->rand,&r));
41c4762a1bSJed Brown   if (r > .9) R = .5;
42c4762a1bSJed Brown   if (r < .1) R = .4;
43c4762a1bSJed Brown   R = .4;
44c4762a1bSJed Brown   /*  The next three lines allow us to access the entries of the vectors directly */
455f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(U,&u));
465f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(Udot,&udot));
475f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(F,&f));
48c4762a1bSJed Brown   f[0] = 2.0*ctx->H*udot[0]/ctx->omega_s + ctx->E*ctx->V*PetscSinScalar(u[1])/ctx->X - R;
49c4762a1bSJed Brown   f[1] = udot[1] - u[0] + ctx->omega_s;
50c4762a1bSJed Brown 
515f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(U,&u));
525f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(Udot,&udot));
535f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(F,&f));
54c4762a1bSJed Brown   PetscFunctionReturn(0);
55c4762a1bSJed Brown }
56c4762a1bSJed Brown 
57c4762a1bSJed Brown /*
58c4762a1bSJed Brown      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
59c4762a1bSJed Brown */
60c4762a1bSJed Brown static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,AppCtx *ctx)
61c4762a1bSJed Brown {
62c4762a1bSJed Brown   PetscInt          rowcol[] = {0,1};
63c4762a1bSJed Brown   PetscScalar       J[2][2];
64c4762a1bSJed Brown   const PetscScalar *u,*udot;
65c4762a1bSJed Brown 
66c4762a1bSJed Brown   PetscFunctionBegin;
675f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(U,&u));
685f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(Udot,&udot));
69c4762a1bSJed Brown   J[0][0] = 2.0*ctx->H*a/ctx->omega_s;   J[0][1] = -ctx->E*ctx->V*PetscCosScalar(u[1])/ctx->X;
70c4762a1bSJed Brown   J[1][0] = -1.0;                        J[1][1] = a;
715f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES));
725f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(U,&u));
735f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(Udot,&udot));
74c4762a1bSJed Brown 
755f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
765f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
77c4762a1bSJed Brown   if (A != B) {
785f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
795f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
80c4762a1bSJed Brown   }
81c4762a1bSJed Brown   PetscFunctionReturn(0);
82c4762a1bSJed Brown }
83c4762a1bSJed Brown 
84c4762a1bSJed Brown int main(int argc,char **argv)
85c4762a1bSJed Brown {
86c4762a1bSJed Brown   TS             ts;            /* ODE integrator */
87c4762a1bSJed Brown   Vec            U;             /* solution will be stored here */
88c4762a1bSJed Brown   Mat            A;             /* Jacobian matrix */
89c4762a1bSJed Brown   PetscErrorCode ierr;
90c4762a1bSJed Brown   PetscMPIInt    size;
91c4762a1bSJed Brown   PetscInt       n = 2;
92c4762a1bSJed Brown   AppCtx         ctx;
93c4762a1bSJed Brown   PetscScalar    *u;
94c4762a1bSJed Brown 
95c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
96c4762a1bSJed Brown      Initialize program
97c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
98*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help));
995f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
1003c633725SBarry Smith   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Only for sequential runs");
101c4762a1bSJed Brown 
102c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
103c4762a1bSJed Brown     Create necessary matrix and vectors
104c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1055f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
1065f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE));
1075f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(A));
1085f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(A));
109c4762a1bSJed Brown 
1105f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(A,&U,NULL));
111c4762a1bSJed Brown 
112c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
113c4762a1bSJed Brown     Set runtime options
114c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
115c4762a1bSJed Brown   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Reaction options","");CHKERRQ(ierr);
116c4762a1bSJed Brown   {
117c4762a1bSJed Brown     ctx.omega_s = 1.0;
1185f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsScalar("-omega_s","","",ctx.omega_s,&ctx.omega_s,NULL));
119c4762a1bSJed Brown     ctx.H       = 1.0;
1205f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsScalar("-H","","",ctx.H,&ctx.H,NULL));
121c4762a1bSJed Brown     ctx.E       = 1.0;
1225f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsScalar("-E","","",ctx.E,&ctx.E,NULL));
123c4762a1bSJed Brown     ctx.V       = 1.0;
1245f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsScalar("-V","","",ctx.V,&ctx.V,NULL));
125c4762a1bSJed Brown     ctx.X       = 1.0;
1265f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsScalar("-X","","",ctx.X,&ctx.X,NULL));
127c4762a1bSJed Brown 
1285f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetArray(U,&u));
129c4762a1bSJed Brown     u[0] = 1;
130c4762a1bSJed Brown     u[1] = .7;
1315f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreArray(U,&u));
1325f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsGetVec(NULL,NULL,"-initial",U,NULL));
133c4762a1bSJed Brown   }
134c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
135c4762a1bSJed Brown 
1365f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&ctx.rand));
1375f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomSetFromOptions(ctx.rand));
138c4762a1bSJed Brown 
139c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
140c4762a1bSJed Brown      Create timestepping solver context
141c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1425f80ce2aSJacob Faibussowitsch   CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts));
1435f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetProblemType(ts,TS_NONLINEAR));
1445f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetType(ts,TSROSW));
1455f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetIFunction(ts,NULL,(TSIFunction) IFunction,&ctx));
1465f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,&ctx));
147c4762a1bSJed Brown 
148c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
149c4762a1bSJed Brown      Set initial conditions
150c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1515f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetSolution(ts,U));
152c4762a1bSJed Brown 
153c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
154c4762a1bSJed Brown      Set solver options
155c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1565f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetMaxTime(ts,2000.0));
1575f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP));
1585f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTimeStep(ts,.001));
1595f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(ts));
160c4762a1bSJed Brown 
161c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
162c4762a1bSJed Brown      Solve nonlinear system
163c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1645f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSolve(ts,U));
165c4762a1bSJed Brown 
166c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
167c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
168c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1695f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
1705f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&U));
1715f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDestroy(&ts));
1725f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomDestroy(&ctx.rand));
173*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
174*b122ec5aSJacob Faibussowitsch   return 0;
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: -ksp_gmres_cgs_refinement_type refine_always -snes_type newtonls -ts_max_steps 10
184c4762a1bSJed Brown 
185c4762a1bSJed Brown    test:
186c4762a1bSJed Brown       suffix: 2
187c4762a1bSJed Brown       args: -ts_max_steps 10
188c4762a1bSJed Brown 
189c4762a1bSJed Brown TEST*/
190