Actual source code: ex3.c

petsc-3.6.4 2016-04-12
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,Pm,E,V,X;
 38:   PetscReal   tf,tcl;
 39: } AppCtx;

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

 53:   /*  The next three lines allow us to access the entries of the vectors directly */
 54:   VecGetArrayRead(U,&u);
 55:   VecGetArrayRead(Udot,&udot);
 56:   VecGetArray(F,&f);
 57:   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 */
 58:   else Pmax = ctx->Pmax;

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

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

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

 82:   VecGetArrayRead(U,&u);
 83:   VecGetArrayRead(Udot,&udot);
 84:   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 */
 85:   else Pmax = ctx->Pmax;

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

 90:   MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);
 91:   VecRestoreArrayRead(U,&u);
 92:   VecRestoreArrayRead(Udot,&udot);

 94:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 95:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 96:   if (A != B) {
 97:     MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 98:     MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
 99:   }
100:   return(0);
101: }

105: int main(int argc,char **argv)
106: {
107:   TS             ts;            /* ODE integrator */
108:   Vec            U;             /* solution will be stored here */
109:   Mat            A;             /* Jacobian matrix */
111:   PetscMPIInt    size;
112:   PetscInt       n = 2;
113:   AppCtx         ctx;
114:   PetscScalar    *u;
115:   PetscReal      du[2] = {0.0,0.0};
116:   PetscBool      ensemble = PETSC_FALSE,flg1,flg2;

118:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
119:      Initialize program
120:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
121:   PetscInitialize(&argc,&argv,(char*)0,help);
122:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
123:   if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs");

125:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
126:     Create necessary matrix and vectors
127:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
128:   MatCreate(PETSC_COMM_WORLD,&A);
129:   MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);
130:   MatSetType(A,MATDENSE);
131:   MatSetFromOptions(A);
132:   MatSetUp(A);

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

136:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
137:     Set runtime options
138:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
139:   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Swing equation options","");
140:   {
141:     ctx.omega_b = 1.0;
142:     ctx.omega_s = 2.0*PETSC_PI*60.0;
143:     ctx.H       = 5.0;
144:     PetscOptionsScalar("-Inertia","","",ctx.H,&ctx.H,NULL);
145:     ctx.D       = 5.0;
146:     PetscOptionsScalar("-D","","",ctx.D,&ctx.D,NULL);
147:     ctx.E       = 1.1378;
148:     ctx.V       = 1.0;
149:     ctx.X       = 0.545;
150:     ctx.Pmax    = ctx.E*ctx.V/ctx.X;;
151:     PetscOptionsScalar("-Pmax","","",ctx.Pmax,&ctx.Pmax,NULL);
152:     ctx.Pm      = 0.9;
153:     PetscOptionsScalar("-Pm","","",ctx.Pm,&ctx.Pm,NULL);
154:     ctx.tf      = 1.0;
155:     ctx.tcl     = 1.05;
156:     PetscOptionsReal("-tf","Time to start fault","",ctx.tf,&ctx.tf,NULL);
157:     PetscOptionsReal("-tcl","Time to end fault","",ctx.tcl,&ctx.tcl,NULL);
158:     PetscOptionsBool("-ensemble","Run ensemble of different initial conditions","",ensemble,&ensemble,NULL);
159:     if (ensemble) {
160:       ctx.tf      = -1;
161:       ctx.tcl     = -1;
162:     }

164:     VecGetArray(U,&u);
165:     u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax);
166:     u[1] = 1.0;
167:     PetscOptionsRealArray("-u","Initial solution","",u,&n,&flg1);
168:     n    = 2;
169:     PetscOptionsRealArray("-du","Perturbation in initial solution","",du,&n,&flg2);
170:     u[0] += du[0];
171:     u[1] += du[1];
172:     VecRestoreArray(U,&u);
173:     if (flg1 || flg2) {
174:       ctx.tf      = -1;
175:       ctx.tcl     = -1;
176:     }
177:   }
178:   PetscOptionsEnd();

180:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
181:      Create timestepping solver context
182:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
183:   TSCreate(PETSC_COMM_WORLD,&ts);
184:   TSSetProblemType(ts,TS_NONLINEAR);
185:   TSSetType(ts,TSTHETA);
186:   TSSetIFunction(ts,NULL,(TSIFunction) IFunction,&ctx);
187:   TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,&ctx);
188:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);

190:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
191:      Set initial conditions
192:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
193:   TSSetSolution(ts,U);

195:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
196:      Set solver options
197:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
198:   TSSetDuration(ts,100000,35.0);
199:   TSSetInitialTimeStep(ts,0.0,.01);
200:   TSSetFromOptions(ts);


203:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
204:      Solve nonlinear system
205:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
206:   if (ensemble) {
207:     for (du[1] = -2.5; du[1] <= .01; du[1] += .1) {
208:       VecGetArray(U,&u);
209:       u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax);
210:       u[1] = ctx.omega_s;
211:       u[0] += du[0];
212:       u[1] += du[1];
213:       VecRestoreArray(U,&u);
214:       TSSetInitialTimeStep(ts,0.0,.01);
215:       TSSolve(ts,U);
216:     }
217:   } else {
218:     TSSolve(ts,U);
219:   }
220:   VecView(U,PETSC_VIEWER_STDOUT_WORLD);
221:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
222:      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
223:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
224:   MatDestroy(&A);
225:   VecDestroy(&U);
226:   TSDestroy(&ts);

228:   PetscFinalize();
229:   return(0);
230: }