Actual source code: ex3.c


  2: static char help[] = "Basic equation for generator stability analysis.\n";

  4: /*F

  6: \begin{eqnarray}
  7:                  \frac{d \theta}{dt} = \omega_b (\omega - \omega_s)
  8:                  \frac{2 H}{\omega_s}\frac{d \omega}{dt} & = & P_m - P_max \sin(\theta) -D(\omega - \omega_s)\\
  9: \end{eqnarray}

 11:   Ensemble of initial conditions
 12:    ./ex3 -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

 14:   Fault at .1 seconds
 15:    ./ex3           -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

 17:   Initial conditions same as when fault is ended
 18:    ./ex3 -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

 20: F*/

 22: /*
 23:    Include "petscts.h" so that we can use TS solvers.  Note that this
 24:    file automatically includes:
 25:      petscsys.h       - base PETSc routines   petscvec.h - vectors
 26:      petscmat.h - matrices
 27:      petscis.h     - index sets            petscksp.h - Krylov subspace methods
 28:      petscviewer.h - viewers               petscpc.h  - preconditioners
 29:      petscksp.h   - linear solvers
 30: */

 32: #include <petscts.h>
 33: #include "ex3.h"

 35: int main(int argc, char **argv)
 36: {
 37:   TS           ts; /* ODE integrator */
 38:   Vec          U;  /* solution will be stored here */
 39:   Mat          A;  /* Jacobian matrix */
 40:   PetscMPIInt  size;
 41:   PetscInt     n = 2;
 42:   AppCtx       ctx;
 43:   PetscScalar *u;
 44:   PetscReal    du[2]    = {0.0, 0.0};
 45:   PetscBool    ensemble = PETSC_FALSE, flg1, flg2;
 46:   PetscInt     direction[2];
 47:   PetscBool    terminate[2];

 49:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 50:      Initialize program
 51:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 52:   PetscFunctionBeginUser;
 53:   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
 54:   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
 55:   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs");

 57:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 58:     Create necessary matrix and vectors
 59:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 60:   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
 61:   PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE));
 62:   PetscCall(MatSetType(A, MATDENSE));
 63:   PetscCall(MatSetFromOptions(A));
 64:   PetscCall(MatSetUp(A));

 66:   PetscCall(MatCreateVecs(A, &U, NULL));

 68:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 69:     Set runtime options
 70:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 71:   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Swing equation options", "");
 72:   {
 73:     ctx.omega_b = 1.0;
 74:     ctx.omega_s = 2.0 * PETSC_PI * 60.0;
 75:     ctx.H       = 5.0;
 76:     PetscCall(PetscOptionsScalar("-Inertia", "", "", ctx.H, &ctx.H, NULL));
 77:     ctx.D = 5.0;
 78:     PetscCall(PetscOptionsScalar("-D", "", "", ctx.D, &ctx.D, NULL));
 79:     ctx.E        = 1.1378;
 80:     ctx.V        = 1.0;
 81:     ctx.X        = 0.545;
 82:     ctx.Pmax     = ctx.E * ctx.V / ctx.X;
 83:     ctx.Pmax_ini = ctx.Pmax;
 84:     PetscCall(PetscOptionsScalar("-Pmax", "", "", ctx.Pmax, &ctx.Pmax, NULL));
 85:     ctx.Pm = 0.9;
 86:     PetscCall(PetscOptionsScalar("-Pm", "", "", ctx.Pm, &ctx.Pm, NULL));
 87:     ctx.tf  = 1.0;
 88:     ctx.tcl = 1.05;
 89:     PetscCall(PetscOptionsReal("-tf", "Time to start fault", "", ctx.tf, &ctx.tf, NULL));
 90:     PetscCall(PetscOptionsReal("-tcl", "Time to end fault", "", ctx.tcl, &ctx.tcl, NULL));
 91:     PetscCall(PetscOptionsBool("-ensemble", "Run ensemble of different initial conditions", "", ensemble, &ensemble, NULL));
 92:     if (ensemble) {
 93:       ctx.tf  = -1;
 94:       ctx.tcl = -1;
 95:     }

 97:     PetscCall(VecGetArray(U, &u));
 98:     u[0] = PetscAsinScalar(ctx.Pm / ctx.Pmax);
 99:     u[1] = 1.0;
100:     PetscCall(PetscOptionsRealArray("-u", "Initial solution", "", u, &n, &flg1));
101:     n = 2;
102:     PetscCall(PetscOptionsRealArray("-du", "Perturbation in initial solution", "", du, &n, &flg2));
103:     u[0] += du[0];
104:     u[1] += du[1];
105:     PetscCall(VecRestoreArray(U, &u));
106:     if (flg1 || flg2) {
107:       ctx.tf  = -1;
108:       ctx.tcl = -1;
109:     }
110:   }
111:   PetscOptionsEnd();

113:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
114:      Create timestepping solver context
115:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
116:   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
117:   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
118:   PetscCall(TSSetType(ts, TSTHETA));
119:   PetscCall(TSSetEquationType(ts, TS_EQ_IMPLICIT));
120:   PetscCall(TSARKIMEXSetFullyImplicit(ts, PETSC_TRUE));
121:   PetscCall(TSSetIFunction(ts, NULL, (TSIFunction)IFunction, &ctx));
122:   PetscCall(TSSetIJacobian(ts, A, A, (TSIJacobian)IJacobian, &ctx));

124:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
125:      Set initial conditions
126:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
127:   PetscCall(TSSetSolution(ts, U));

129:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
130:      Set solver options
131:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
132:   PetscCall(TSSetMaxTime(ts, 35.0));
133:   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
134:   PetscCall(TSSetTimeStep(ts, .1));
135:   PetscCall(TSSetFromOptions(ts));

137:   direction[0] = direction[1] = 1;
138:   terminate[0] = terminate[1] = PETSC_FALSE;

140:   PetscCall(TSSetEventHandler(ts, 2, direction, terminate, EventFunction, PostEventFunction, (void *)&ctx));

142:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
143:      Solve nonlinear system
144:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
145:   if (ensemble) {
146:     for (du[1] = -2.5; du[1] <= .01; du[1] += .1) {
147:       PetscCall(VecGetArray(U, &u));
148:       u[0] = PetscAsinScalar(ctx.Pm / ctx.Pmax);
149:       u[1] = ctx.omega_s;
150:       u[0] += du[0];
151:       u[1] += du[1];
152:       PetscCall(VecRestoreArray(U, &u));
153:       PetscCall(TSSetTimeStep(ts, .01));
154:       PetscCall(TSSolve(ts, U));
155:     }
156:   } else {
157:     PetscCall(TSSolve(ts, U));
158:   }
159:   PetscCall(VecView(U, PETSC_VIEWER_STDOUT_WORLD));
160:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
161:      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
162:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
163:   PetscCall(MatDestroy(&A));
164:   PetscCall(VecDestroy(&U));
165:   PetscCall(TSDestroy(&ts));
166:   PetscCall(PetscFinalize());
167:   return 0;
168: }

170: /*TEST

172:    build:
173:      requires: !complex !single

175:    test:
176:       args: -nox

178: TEST*/