Actual source code: ex3.c


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


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

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

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

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

 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 */
 41:   PetscMPIInt    size;
 42:   PetscInt       n = 2;
 43:   AppCtx         ctx;
 44:   PetscScalar    *u;
 45:   PetscReal      du[2] = {0.0,0.0};
 46:   PetscBool      ensemble = PETSC_FALSE,flg1,flg2;
 47:   PetscInt       direction[2];
 48:   PetscBool      terminate[2];

 50:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 51:      Initialize program
 52:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 53:   PetscInitialize(&argc,&argv,(char*)0,help);
 54:   MPI_Comm_size(PETSC_COMM_WORLD,&size);

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

 66:   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:     PetscOptionsScalar("-Inertia","","",ctx.H,&ctx.H,NULL);
 77:     ctx.D       = 5.0;
 78:     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:     PetscOptionsScalar("-Pmax","","",ctx.Pmax,&ctx.Pmax,NULL);
 85:     ctx.Pm      = 0.9;
 86:     PetscOptionsScalar("-Pm","","",ctx.Pm,&ctx.Pm,NULL);
 87:     ctx.tf      = 1.0;
 88:     ctx.tcl     = 1.05;
 89:     PetscOptionsReal("-tf","Time to start fault","",ctx.tf,&ctx.tf,NULL);
 90:     PetscOptionsReal("-tcl","Time to end fault","",ctx.tcl,&ctx.tcl,NULL);
 91:     PetscOptionsBool("-ensemble","Run ensemble of different initial conditions","",ensemble,&ensemble,NULL);
 92:     if (ensemble) {
 93:       ctx.tf      = -1;
 94:       ctx.tcl     = -1;
 95:     }

 97:     VecGetArray(U,&u);
 98:     u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax);
 99:     u[1] = 1.0;
100:     PetscOptionsRealArray("-u","Initial solution","",u,&n,&flg1);
101:     n    = 2;
102:     PetscOptionsRealArray("-du","Perturbation in initial solution","",du,&n,&flg2);
103:     u[0] += du[0];
104:     u[1] += du[1];
105:     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:   TSCreate(PETSC_COMM_WORLD,&ts);
117:   TSSetProblemType(ts,TS_NONLINEAR);
118:   TSSetType(ts,TSTHETA);
119:   TSSetEquationType(ts,TS_EQ_IMPLICIT);
120:   TSARKIMEXSetFullyImplicit(ts,PETSC_TRUE);
121:   TSSetIFunction(ts,NULL,(TSIFunction) IFunction,&ctx);
122:   TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,&ctx);

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

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

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

140:   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:       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:       VecRestoreArray(U,&u);
153:       TSSetTimeStep(ts,.01);
154:       TSSolve(ts,U);
155:     }
156:   } else {
157:     TSSolve(ts,U);
158:   }
159:   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:   MatDestroy(&A);
164:   VecDestroy(&U);
165:   TSDestroy(&ts);
166:   PetscFinalize();
167:   return 0;
168: }

170: /*TEST

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

175:    test:
176:       args: -nox

178: TEST*/