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: */
 31: /*T

 33: T*/

 35: #include <petscts.h>
 36: #include "ex3.h"

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

 53:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 54:      Initialize program
 55:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 56:   PetscInitialize(&argc,&argv,(char*)0,help);
 57:   MPI_Comm_size(PETSC_COMM_WORLD,&size);

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

 69:   MatCreateVecs(A,&U,NULL);

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

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

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

127:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
128:      Set initial conditions
129:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
130:   TSSetSolution(ts,U);

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

140:   direction[0] = direction[1] = 1;
141:   terminate[0] = terminate[1] = PETSC_FALSE;

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

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

173: /*TEST

175:    build:
176:      requires: !complex !single

178:    test:
179:       args: -nox

181: TEST*/