Actual source code: ex2.c

petsc-3.8.3 2017-12-09
Report Typos and Errors

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


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



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


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

 36: typedef struct {
 37:   PetscScalar H,D,omega_s,Pmax,Pm,E,V,X;
 38:   PetscReal   tf,tcl;
 39: } AppCtx;

 41: /*
 42:      Defines the ODE passed to the ODE solver
 43: */
 44: static PetscErrorCode IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,AppCtx *ctx)
 45: {
 46:   PetscErrorCode    ierr;
 47:   PetscScalar       *f,Pmax;
 48:   const PetscScalar *u,*udot;

 51:   /*  The next three lines allow us to access the entries of the vectors directly */
 52:   VecGetArrayRead(U,&u);
 53:   VecGetArrayRead(Udot,&udot);
 54:   VecGetArray(F,&f);
 55:   if ((t > ctx->tf) && (t < ctx->tcl)) Pmax = 0.0; /* A short-circuit on the generator terminal that drives the electrical power output (Pmax*sin(delta)) to 0 */
 56:   else if (t >= ctx->tcl) Pmax = ctx->E/0.745;
 57:   else Pmax = ctx->Pmax;
 58:   f[0] = udot[0] - ctx->omega_s*(u[1] - 1.0);
 59:   f[1] = 2.0*ctx->H*udot[1] +  Pmax*PetscSinScalar(u[0]) + ctx->D*(u[1] - 1.0)- ctx->Pm;

 61:   VecRestoreArrayRead(U,&u);
 62:   VecRestoreArrayRead(Udot,&udot);
 63:   VecRestoreArray(F,&f);
 64:   return(0);
 65: }

 67: /*
 68:      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
 69: */
 70: static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,AppCtx *ctx)
 71: {
 72:   PetscErrorCode    ierr;
 73:   PetscInt          rowcol[] = {0,1};
 74:   PetscScalar       J[2][2],Pmax;
 75:   const PetscScalar *u,*udot;

 78:   VecGetArrayRead(U,&u);
 79:   VecGetArrayRead(Udot,&udot);
 80:   if ((t > ctx->tf) && (t < ctx->tcl)) Pmax = 0.0; /* A short-circuit on the generator terminal that drives the electrical power output (Pmax*sin(delta)) to 0 */
 81:   else if (t >= ctx->tcl) Pmax = ctx->E/0.745;
 82:   else Pmax = ctx->Pmax;

 84:   J[0][0] = a;                       J[0][1] = -ctx->omega_s;
 85:   J[1][1] = 2.0*ctx->H*a + ctx->D;   J[1][0] = Pmax*PetscCosScalar(u[0]);

 87:   MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);
 88:   VecRestoreArrayRead(U,&u);
 89:   VecRestoreArrayRead(Udot,&udot);

 91:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 92:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 93:   if (A != B) {
 94:     MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 95:     MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
 96:   }
 97:   return(0);
 98: }


101: PetscErrorCode PostStep(TS ts)
102: {
104:   Vec            X;
105:   PetscReal      t;

108:   TSGetTime(ts,&t);
109:   if (t >= .2) {
110:     TSGetSolution(ts,&X);
111:     VecView(X,PETSC_VIEWER_STDOUT_WORLD);
112:     exit(0);
113:     /* results in initial conditions after fault of -u 0.496792,1.00932 */
114:   }
115:   return(0);
116: }


119: int main(int argc,char **argv)
120: {
121:   TS             ts;            /* ODE integrator */
122:   Vec            U;             /* solution will be stored here */
123:   Mat            A;             /* Jacobian matrix */
125:   PetscMPIInt    size;
126:   PetscInt       n = 2;
127:   AppCtx         ctx;
128:   PetscScalar    *u;
129:   PetscReal      du[2] = {0.0,0.0};
130:   PetscBool      ensemble = PETSC_FALSE,flg1,flg2;

132:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
133:      Initialize program
134:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
135:   PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
136:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
137:   if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs");

139:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
140:     Create necessary matrix and vectors
141:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
142:   MatCreate(PETSC_COMM_WORLD,&A);
143:   MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);
144:   MatSetType(A,MATDENSE);
145:   MatSetFromOptions(A);
146:   MatSetUp(A);

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

150:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
151:     Set runtime options
152:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
153:   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Swing equation options","");
154:   {
155:     ctx.omega_s = 2.0*PETSC_PI*60.0;
156:     ctx.H       = 5.0;
157:     PetscOptionsScalar("-Inertia","","",ctx.H,&ctx.H,NULL);
158:     ctx.D       = 5.0;
159:     PetscOptionsScalar("-D","","",ctx.D,&ctx.D,NULL);
160:     ctx.E       = 1.1378;
161:     ctx.V       = 1.0;
162:     ctx.X       = 0.545;
163:     ctx.Pmax    = ctx.E*ctx.V/ctx.X;;
164:     PetscOptionsScalar("-Pmax","","",ctx.Pmax,&ctx.Pmax,NULL);
165:     ctx.Pm      = 0.9;
166:     PetscOptionsScalar("-Pm","","",ctx.Pm,&ctx.Pm,NULL);
167:     ctx.tf      = 1.0;
168:     ctx.tcl     = 1.05;
169:     PetscOptionsReal("-tf","Time to start fault","",ctx.tf,&ctx.tf,NULL);
170:     PetscOptionsReal("-tcl","Time to end fault","",ctx.tcl,&ctx.tcl,NULL);
171:     PetscOptionsBool("-ensemble","Run ensemble of different initial conditions","",ensemble,&ensemble,NULL);
172:     if (ensemble) {
173:       ctx.tf      = -1;
174:       ctx.tcl     = -1;
175:     }

177:     VecGetArray(U,&u);
178:     u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax);
179:     u[1] = 1.0;
180:     PetscOptionsRealArray("-u","Initial solution","",u,&n,&flg1);
181:     n    = 2;
182:     PetscOptionsRealArray("-du","Perturbation in initial solution","",du,&n,&flg2);
183:     u[0] += du[0];
184:     u[1] += du[1];
185:     VecRestoreArray(U,&u);
186:     if (flg1 || flg2) {
187:       ctx.tf      = -1;
188:       ctx.tcl     = -1;
189:     }
190:   }
191:   PetscOptionsEnd();

193:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
194:      Create timestepping solver context
195:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
196:   TSCreate(PETSC_COMM_WORLD,&ts);
197:   TSSetProblemType(ts,TS_NONLINEAR);
198:   TSSetType(ts,TSROSW);
199:   TSSetIFunction(ts,NULL,(TSIFunction) IFunction,&ctx);
200:   TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,&ctx);

202:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
203:      Set initial conditions
204:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
205:   TSSetSolution(ts,U);

207:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
208:      Set solver options
209:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
210:   TSSetMaxTime(ts,35.0);
211:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);
212:   TSSetTimeStep(ts,.01);
213:   TSSetFromOptions(ts);
214:   /* TSSetPostStep(ts,PostStep);  */


217:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
218:      Solve nonlinear system
219:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
220:   if (ensemble) {
221:     for (du[1] = -2.5; du[1] <= .01; du[1] += .1) {
222:       VecGetArray(U,&u);
223:       u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax);
224:       u[1] = ctx.omega_s;
225:       u[0] += du[0];
226:       u[1] += du[1];
227:       VecRestoreArray(U,&u);
228:       TSSetTimeStep(ts,.01);
229:       TSSolve(ts,U);
230:     }
231:   } else {
232:     TSSolve(ts,U);
233:   }
234:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
235:      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
236:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
237:   MatDestroy(&A);
238:   VecDestroy(&U);
239:   TSDestroy(&ts);
240:   PetscFinalize();
241:   return ierr;
242: }