Actual source code: ex2.c petsc-3.6.2 2015-10-02 
   
  2:  static char help[] = "Basic equation for generator stability analysis.\n" ;
 
                 
 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;
 43:  /* 
 44:       Defines the ODE passed to the ODE solver 
 45:  */ 
 46:  PetscErrorCode  IFunction(TS  ts,PetscReal  t,Vec  U,Vec  Udot,Vec  F,AppCtx *ctx) 47:  {
 48:    PetscErrorCode     ierr;
 49:    PetscScalar        *f,Pmax;
 50:    const PetscScalar  *u,*udot;
 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  if  (t >= ctx->tcl) Pmax = ctx->E/0.745;
 59:    else  Pmax = ctx->Pmax;
 60:    f[0] = udot[0] - ctx->omega_s*(u[1] - 1.0);
 61:    f[1] = 2.0*ctx->H*udot[1] +  Pmax*PetscSinScalar(u[0]) + ctx->D*(u[1] - 1.0)- 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:  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  if  (t >= ctx->tcl) Pmax = ctx->E/0.745;
 86:    else  Pmax = ctx->Pmax;
 88:    J[0][0] = a;                       J[0][1] = -ctx->omega_s;
 89:    J[1][1] = 2.0*ctx->H*a + ctx->D;   J[1][0] = Pmax*PetscCosScalar(u[0]);
 91:    MatSetValues (B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES );
 92:    VecRestoreArrayRead (U,&u);
 93:    VecRestoreArrayRead (Udot,&udot);
 95:    MatAssemblyBegin (A,MAT_FINAL_ASSEMBLY);
 96:    MatAssemblyEnd (A,MAT_FINAL_ASSEMBLY);
 97:    if  (A != B) {
 98:      MatAssemblyBegin (B,MAT_FINAL_ASSEMBLY);
 99:      MatAssemblyEnd (B,MAT_FINAL_ASSEMBLY);
100:    }
101:    return (0);
102:  }
107:  PetscErrorCode  PostStep(TS  ts)108:  {
110:    Vec             X;
111:    PetscReal       t;
114:    TSGetTime (ts,&t);
115:    if  (t >= .2) {
116:      TSGetSolution (ts,&X);
117:      VecView (X,PETSC_VIEWER_STDOUT_WORLD );
118:      exit(0);
119:      /* results in initial conditions after fault of -u 0.496792,1.00932 */ 
120:    }
121:    return (0);
122:  }
127:  128:  {
129:    TS              ts;            /* ODE integrator */ 
130:    Vec             U;             /* solution will be stored here */ 
131:    Mat             A;             /* Jacobian matrix */ 
133:    PetscMPIInt     size;
134:    PetscInt        n = 2;
135:    AppCtx         ctx;
136:    PetscScalar     *u;
137:    PetscReal       du[2] = {0.0,0.0};
138:    PetscBool       ensemble = PETSC_FALSE ,flg1,flg2;
140:    /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
141:       Initialize program 
142:       - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 
143:    PetscInitialize (&argc,&argv,(char*)0,help);
144:    MPI_Comm_size (PETSC_COMM_WORLD ,&size);
145:    if  (size > 1) SETERRQ (PETSC_COMM_WORLD ,PETSC_ERR_SUP,"Only for sequential runs" );
147:    /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
148:      Create necessary matrix and vectors 
149:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 
150:    MatCreate (PETSC_COMM_WORLD ,&A);
151:    MatSetSizes (A,n,n,PETSC_DETERMINE ,PETSC_DETERMINE );
152:    MatSetType (A,MATDENSE );
153:    MatSetFromOptions (A);
154:    MatSetUp (A);
156:    MatCreateVecs (A,&U,NULL);
158:    /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
159:      Set runtime options 
160:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 
161:    PetscOptionsBegin (PETSC_COMM_WORLD ,NULL,"Swing equation options" ,"" );
162:    {
163:      ctx.omega_s = 2.0*PETSC_PI*60.0;
164:      ctx.H       = 5.0;
165:      PetscOptionsScalar ("-Inertia" ,"" ,"" ,ctx.H,&ctx.H,NULL);
166:      ctx.D       = 5.0;
167:      PetscOptionsScalar ("-D" ,"" ,"" ,ctx.D,&ctx.D,NULL);
168:      ctx.E       = 1.1378;
169:      ctx.V       = 1.0;
170:      ctx.X       = 0.545;
171:      ctx.Pmax    = ctx.E*ctx.V/ctx.X;;
172:      PetscOptionsScalar ("-Pmax" ,"" ,"" ,ctx.Pmax,&ctx.Pmax,NULL);
173:      ctx.Pm      = 0.9;
174:      PetscOptionsScalar ("-Pm" ,"" ,"" ,ctx.Pm,&ctx.Pm,NULL);
175:      ctx.tf      = 1.0;
176:      ctx.tcl     = 1.05;
177:      PetscOptionsReal ("-tf" ,"Time to start fault" ,"" ,ctx.tf,&ctx.tf,NULL);
178:      PetscOptionsReal ("-tcl" ,"Time to end fault" ,"" ,ctx.tcl,&ctx.tcl,NULL);
179:      PetscOptionsBool ("-ensemble" ,"Run ensemble of different initial conditions" ,"" ,ensemble,&ensemble,NULL);
180:      if  (ensemble) {
181:        ctx.tf      = -1;
182:        ctx.tcl     = -1;
183:      }
185:      VecGetArray (U,&u);
186:      u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax);
187:      u[1] = 1.0;
188:      PetscOptionsRealArray ("-u" ,"Initial solution" ,"" ,u,&n,&flg1);
189:      n    = 2;
190:      PetscOptionsRealArray ("-du" ,"Perturbation in initial solution" ,"" ,du,&n,&flg2);
191:      u[0] += du[0];
192:      u[1] += du[1];
193:      VecRestoreArray (U,&u);
194:      if  (flg1 || flg2) {
195:        ctx.tf      = -1;
196:        ctx.tcl     = -1;
197:      }
198:    }
199:    PetscOptionsEnd ();
201:    /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
202:       Create timestepping solver context 
203:       - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 
204:    TSCreate (PETSC_COMM_WORLD ,&ts);
205:    TSSetProblemType (ts,TS_NONLINEAR);
206:    TSSetType (ts,TSROSW );
207:    TSSetIFunction (ts,NULL,(TSIFunction) IFunction,&ctx);
208:    TSSetIJacobian (ts,A,A,(TSIJacobian)IJacobian,&ctx);
210:    /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
211:       Set initial conditions 
212:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 
213:    TSSetSolution (ts,U);
215:    /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
216:       Set solver options 
217:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 
218:    TSSetDuration (ts,100000,35.0);
219:    TSSetInitialTimeStep (ts,0.0,.01);
220:    TSSetFromOptions (ts);
221:    /* TSSetPostStep (ts,PostStep);  */ 
224:    /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
225:       Solve nonlinear system 
226:       - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 
227:    if  (ensemble) {
228:      for  (du[1] = -2.5; du[1] <= .01; du[1] += .1) {
229:        VecGetArray (U,&u);
230:        u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax);
231:        u[1] = ctx.omega_s;
232:        u[0] += du[0];
233:        u[1] += du[1];
234:        VecRestoreArray (U,&u);
235:        TSSetInitialTimeStep (ts,0.0,.01);
236:        TSSolve (ts,U);
237:      }
238:    } else  {
239:      TSSolve (ts,U);
240:    }
241:    /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
242:       Free work space.  All PETSc objects should be destroyed when they are no longer needed. 
243:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 
244:    MatDestroy (&A);
245:    VecDestroy (&U);
246:    TSDestroy (&ts);
248:    PetscFinalize ();
249:    return (0);
250:  }