Actual source code: ex3.c

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

  3: /*F

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

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

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

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

 19: F*/

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

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

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

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

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

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

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

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

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

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

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

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

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

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

169: /*TEST

171:    build:
172:      requires: !complex !single

174:    test:
175:       args: -nox

177: TEST*/