Actual source code: ex3.c

petsc-3.7.4 2016-10-02
Report Typos and Errors
  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
./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_b,omega_s,Pmax,Pmax_ini,Pm,E,V,X;
 38:   PetscReal   tf,tcl;
 39: } AppCtx;

 41: /* Event check */
 44: PetscErrorCode EventFunction(TS ts,PetscReal t,Vec X,PetscScalar *fvalue,void *ctx)
 45: {
 46:   AppCtx        *user=(AppCtx*)ctx;

 49:   /* Event for fault-on time */
 50:   fvalue[0] = t - user->tf;
 51:   /* Event for fault-off time */
 52:   fvalue[1] = t - user->tcl;

 54:   return(0);
 55: }

 59: PetscErrorCode PostEventFunction(TS ts,PetscInt nevents,PetscInt event_list[],PetscReal t,Vec X,PetscBool forwardsolve,void* ctx)
 60: {
 61:   AppCtx *user=(AppCtx*)ctx;
 62: 

 65:   if (event_list[0] == 0) user->Pmax = 0.0; /* Apply disturbance - this is done by setting Pmax = 0 */
 66:   else if(event_list[0] == 1) user->Pmax = user->Pmax_ini; /* Remove the fault  - this is done by setting Pmax = Pmax_ini */
 67:   return(0);
 68: }

 72: /*
 73:      Defines the ODE passed to the ODE solver
 74: */
 75: static PetscErrorCode IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,AppCtx *ctx)
 76: {
 77:   PetscErrorCode    ierr;
 78:   const PetscScalar *u,*udot;
 79:   PetscScalar       *f,Pmax;

 82:   /*  The next three lines allow us to access the entries of the vectors directly */
 83:   VecGetArrayRead(U,&u);
 84:   VecGetArrayRead(Udot,&udot);
 85:   VecGetArray(F,&f);
 86:   Pmax = ctx->Pmax;

 88:   f[0] = udot[0] - ctx->omega_b*(u[1] - ctx->omega_s);
 89:   f[1] = 2.0*ctx->H/ctx->omega_s*udot[1] +  Pmax*PetscSinScalar(u[0]) + ctx->D*(u[1] - ctx->omega_s)- ctx->Pm;

 91:   VecRestoreArrayRead(U,&u);
 92:   VecRestoreArrayRead(Udot,&udot);
 93:   VecRestoreArray(F,&f);
 94:   return(0);
 95: }

 99: /*
100:      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
101: */
102: static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,AppCtx *ctx)
103: {
104:   PetscErrorCode    ierr;
105:   PetscInt          rowcol[] = {0,1};
106:   PetscScalar       J[2][2],Pmax;
107:   const PetscScalar *u,*udot;

110:   VecGetArrayRead(U,&u);
111:   VecGetArrayRead(Udot,&udot);
112:   Pmax = ctx->Pmax;

114:   J[0][0] = a;                                    J[0][1] = -ctx->omega_b;
115:   J[1][1] = 2.0*ctx->H/ctx->omega_s*a + ctx->D;   J[1][0] = Pmax*PetscCosScalar(u[0]);

117:   MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);
118:   VecRestoreArrayRead(U,&u);
119:   VecRestoreArrayRead(Udot,&udot);

121:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
122:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
123:   if (A != B) {
124:     MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
125:     MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
126:   }
127:   return(0);
128: }

132: int main(int argc,char **argv)
133: {
134:   TS             ts;            /* ODE integrator */
135:   Vec            U;             /* solution will be stored here */
136:   Mat            A;             /* Jacobian matrix */
138:   PetscMPIInt    size;
139:   PetscInt       n = 2;
140:   AppCtx         ctx;
141:   PetscScalar    *u;
142:   PetscReal      du[2] = {0.0,0.0};
143:   PetscBool      ensemble = PETSC_FALSE,flg1,flg2;
144:   PetscInt       direction[2];
145:   PetscBool      terminate[2];

147:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
148:      Initialize program
149:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
150:   PetscInitialize(&argc,&argv,(char*)0,help);
151:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
152:   if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs");

154:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
155:     Create necessary matrix and vectors
156:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
157:   MatCreate(PETSC_COMM_WORLD,&A);
158:   MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);
159:   MatSetType(A,MATDENSE);
160:   MatSetFromOptions(A);
161:   MatSetUp(A);

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

165:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
166:     Set runtime options
167:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
168:   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Swing equation options","");
169:   {
170:     ctx.omega_b = 1.0;
171:     ctx.omega_s = 2.0*PETSC_PI*60.0;
172:     ctx.H       = 5.0;
173:     PetscOptionsScalar("-Inertia","","",ctx.H,&ctx.H,NULL);
174:     ctx.D       = 5.0;
175:     PetscOptionsScalar("-D","","",ctx.D,&ctx.D,NULL);
176:     ctx.E       = 1.1378;
177:     ctx.V       = 1.0;
178:     ctx.X       = 0.545;
179:     ctx.Pmax    = ctx.E*ctx.V/ctx.X;
180:     ctx.Pmax_ini = ctx.Pmax;
181:     PetscOptionsScalar("-Pmax","","",ctx.Pmax,&ctx.Pmax,NULL);
182:     ctx.Pm      = 0.9;
183:     PetscOptionsScalar("-Pm","","",ctx.Pm,&ctx.Pm,NULL);
184:     ctx.tf      = 1.0;
185:     ctx.tcl     = 1.05;
186:     PetscOptionsReal("-tf","Time to start fault","",ctx.tf,&ctx.tf,NULL);
187:     PetscOptionsReal("-tcl","Time to end fault","",ctx.tcl,&ctx.tcl,NULL);
188:     PetscOptionsBool("-ensemble","Run ensemble of different initial conditions","",ensemble,&ensemble,NULL);
189:     if (ensemble) {
190:       ctx.tf      = -1;
191:       ctx.tcl     = -1;
192:     }

194:     VecGetArray(U,&u);
195:     u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax);
196:     u[1] = 1.0;
197:     PetscOptionsRealArray("-u","Initial solution","",u,&n,&flg1);
198:     n    = 2;
199:     PetscOptionsRealArray("-du","Perturbation in initial solution","",du,&n,&flg2);
200:     u[0] += du[0];
201:     u[1] += du[1];
202:     VecRestoreArray(U,&u);
203:     if (flg1 || flg2) {
204:       ctx.tf      = -1;
205:       ctx.tcl     = -1;
206:     }
207:   }
208:   PetscOptionsEnd();

210:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
211:      Create timestepping solver context
212:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
213:   TSCreate(PETSC_COMM_WORLD,&ts);
214:   TSSetProblemType(ts,TS_NONLINEAR);
215:   TSSetType(ts,TSTHETA);
216:   TSSetEquationType(ts,TS_EQ_IMPLICIT);
217:   TSARKIMEXSetFullyImplicit(ts,PETSC_TRUE);
218:   TSSetIFunction(ts,NULL,(TSIFunction) IFunction,&ctx);
219:   TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,&ctx);
220:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);

222:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
223:      Set initial conditions
224:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
225:   TSSetSolution(ts,U);

227:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
228:      Set solver options
229:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
230:   TSSetDuration(ts,100000,35.0);
231:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
232:   TSSetInitialTimeStep(ts,0.0,.01);
233:   TSSetFromOptions(ts);

235:   direction[0] = direction[1] = 1;
236:   terminate[0] = terminate[1] = PETSC_FALSE;

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

240:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
241:      Solve nonlinear system
242:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
243:   if (ensemble) {
244:     for (du[1] = -2.5; du[1] <= .01; du[1] += .1) {
245:       VecGetArray(U,&u);
246:       u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax);
247:       u[1] = ctx.omega_s;
248:       u[0] += du[0];
249:       u[1] += du[1];
250:       VecRestoreArray(U,&u);
251:       TSSetInitialTimeStep(ts,0.0,.01);
252:       TSSolve(ts,U);
253:     }
254:   } else {
255:     TSSolve(ts,U);
256:   }
257:   VecView(U,PETSC_VIEWER_STDOUT_WORLD);
258:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
259:      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
260:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
261:   MatDestroy(&A);
262:   VecDestroy(&U);
263:   TSDestroy(&ts);

265:   PetscFinalize();
266:   return(0);
267: }