Actual source code: ex3.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{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 */
 42: PetscErrorCode EventFunction(TS ts,PetscReal t,Vec X,PetscScalar *fvalue,void *ctx)
 43: {
 44:   AppCtx        *user=(AppCtx*)ctx;

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

 52:   return(0);
 53: }

 55: PetscErrorCode PostEventFunction(TS ts,PetscInt nevents,PetscInt event_list[],PetscReal t,Vec X,PetscBool forwardsolve,void* ctx)
 56: {
 57:   AppCtx *user=(AppCtx*)ctx;
 58: 

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

 66: /*
 67:      Defines the ODE passed to the ODE solver
 68: */
 69: static PetscErrorCode IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,AppCtx *ctx)
 70: {
 71:   PetscErrorCode    ierr;
 72:   const PetscScalar *u,*udot;
 73:   PetscScalar       *f,Pmax;

 76:   /*  The next three lines allow us to access the entries of the vectors directly */
 77:   VecGetArrayRead(U,&u);
 78:   VecGetArrayRead(Udot,&udot);
 79:   VecGetArray(F,&f);
 80:   Pmax = ctx->Pmax;

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

 85:   VecRestoreArrayRead(U,&u);
 86:   VecRestoreArrayRead(Udot,&udot);
 87:   VecRestoreArray(F,&f);
 88:   return(0);
 89: }

 91: /*
 92:      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
 93: */
 94: static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,AppCtx *ctx)
 95: {
 96:   PetscErrorCode    ierr;
 97:   PetscInt          rowcol[] = {0,1};
 98:   PetscScalar       J[2][2],Pmax;
 99:   const PetscScalar *u,*udot;

102:   VecGetArrayRead(U,&u);
103:   VecGetArrayRead(Udot,&udot);
104:   Pmax = ctx->Pmax;

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

109:   MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);
110:   VecRestoreArrayRead(U,&u);
111:   VecRestoreArrayRead(Udot,&udot);

113:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
114:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
115:   if (A != B) {
116:     MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
117:     MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
118:   }
119:   return(0);
120: }

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

137:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
138:      Initialize program
139:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
140:   PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
141:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
142:   if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs");

144:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
145:     Create necessary matrix and vectors
146:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
147:   MatCreate(PETSC_COMM_WORLD,&A);
148:   MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);
149:   MatSetType(A,MATDENSE);
150:   MatSetFromOptions(A);
151:   MatSetUp(A);

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

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

184:     VecGetArray(U,&u);
185:     u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax);
186:     u[1] = 1.0;
187:     PetscOptionsRealArray("-u","Initial solution","",u,&n,&flg1);
188:     n    = 2;
189:     PetscOptionsRealArray("-du","Perturbation in initial solution","",du,&n,&flg2);
190:     u[0] += du[0];
191:     u[1] += du[1];
192:     VecRestoreArray(U,&u);
193:     if (flg1 || flg2) {
194:       ctx.tf      = -1;
195:       ctx.tcl     = -1;
196:     }
197:   }
198:   PetscOptionsEnd();

200:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
201:      Create timestepping solver context
202:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
203:   TSCreate(PETSC_COMM_WORLD,&ts);
204:   TSSetProblemType(ts,TS_NONLINEAR);
205:   TSSetType(ts,TSTHETA);
206:   TSSetEquationType(ts,TS_EQ_IMPLICIT);
207:   TSARKIMEXSetFullyImplicit(ts,PETSC_TRUE);
208:   TSSetIFunction(ts,NULL,(TSIFunction) IFunction,&ctx);
209:   TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,&ctx);

211:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
212:      Set initial conditions
213:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
214:   TSSetSolution(ts,U);

216:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
217:      Set solver options
218:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
219:   TSSetMaxTime(ts,35.0);
220:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);
221:   TSSetTimeStep(ts,.1);
222:   TSSetFromOptions(ts);

224:   direction[0] = direction[1] = 1;
225:   terminate[0] = terminate[1] = PETSC_FALSE;

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

229:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
230:      Solve nonlinear system
231:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
232:   if (ensemble) {
233:     for (du[1] = -2.5; du[1] <= .01; du[1] += .1) {
234:       VecGetArray(U,&u);
235:       u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax);
236:       u[1] = ctx.omega_s;
237:       u[0] += du[0];
238:       u[1] += du[1];
239:       VecRestoreArray(U,&u);
240:       TSSetTimeStep(ts,.01);
241:       TSSolve(ts,U);
242:     }
243:   } else {
244:     TSSolve(ts,U);
245:   }
246:   VecView(U,PETSC_VIEWER_STDOUT_WORLD);
247:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
248:      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
249:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
250:   MatDestroy(&A);
251:   VecDestroy(&U);
252:   TSDestroy(&ts);
253:   PetscFinalize();
254:   return ierr;
255: }