Actual source code: ex3.c

petsc-3.12.1 2019-10-22
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
./ex3 -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
./ex3 -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
./ex3 -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: /*T

 36: T*/

 38: #include <petscts.h>
 39: #include "ex3.h"

 41: int main(int argc,char **argv)
 42: {
 43:   TS             ts;            /* ODE integrator */
 44:   Vec            U;             /* solution will be stored here */
 45:   Mat            A;             /* Jacobian matrix */
 47:   PetscMPIInt    size;
 48:   PetscInt       n = 2;
 49:   AppCtx         ctx;
 50:   PetscScalar    *u;
 51:   PetscReal      du[2] = {0.0,0.0};
 52:   PetscBool      ensemble = PETSC_FALSE,flg1,flg2;
 53:   PetscInt       direction[2];
 54:   PetscBool      terminate[2];

 56:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 57:      Initialize program
 58:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 59:   PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 60:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 61:   if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs");

 63:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 64:     Create necessary matrix and vectors
 65:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 66:   MatCreate(PETSC_COMM_WORLD,&A);
 67:   MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);
 68:   MatSetType(A,MATDENSE);
 69:   MatSetFromOptions(A);
 70:   MatSetUp(A);

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

 74:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 75:     Set runtime options
 76:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 77:   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Swing equation options","");
 78:   {
 79:     ctx.omega_b = 1.0;
 80:     ctx.omega_s = 2.0*PETSC_PI*60.0;
 81:     ctx.H       = 5.0;
 82:     PetscOptionsScalar("-Inertia","","",ctx.H,&ctx.H,NULL);
 83:     ctx.D       = 5.0;
 84:     PetscOptionsScalar("-D","","",ctx.D,&ctx.D,NULL);
 85:     ctx.E       = 1.1378;
 86:     ctx.V       = 1.0;
 87:     ctx.X       = 0.545;
 88:     ctx.Pmax    = ctx.E*ctx.V/ctx.X;
 89:     ctx.Pmax_ini = ctx.Pmax;
 90:     PetscOptionsScalar("-Pmax","","",ctx.Pmax,&ctx.Pmax,NULL);
 91:     ctx.Pm      = 0.9;
 92:     PetscOptionsScalar("-Pm","","",ctx.Pm,&ctx.Pm,NULL);
 93:     ctx.tf      = 1.0;
 94:     ctx.tcl     = 1.05;
 95:     PetscOptionsReal("-tf","Time to start fault","",ctx.tf,&ctx.tf,NULL);
 96:     PetscOptionsReal("-tcl","Time to end fault","",ctx.tcl,&ctx.tcl,NULL);
 97:     PetscOptionsBool("-ensemble","Run ensemble of different initial conditions","",ensemble,&ensemble,NULL);
 98:     if (ensemble) {
 99:       ctx.tf      = -1;
100:       ctx.tcl     = -1;
101:     }

103:     VecGetArray(U,&u);
104:     u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax);
105:     u[1] = 1.0;
106:     PetscOptionsRealArray("-u","Initial solution","",u,&n,&flg1);
107:     n    = 2;
108:     PetscOptionsRealArray("-du","Perturbation in initial solution","",du,&n,&flg2);
109:     u[0] += du[0];
110:     u[1] += du[1];
111:     VecRestoreArray(U,&u);
112:     if (flg1 || flg2) {
113:       ctx.tf      = -1;
114:       ctx.tcl     = -1;
115:     }
116:   }
117:   PetscOptionsEnd();

119:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
120:      Create timestepping solver context
121:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
122:   TSCreate(PETSC_COMM_WORLD,&ts);
123:   TSSetProblemType(ts,TS_NONLINEAR);
124:   TSSetType(ts,TSTHETA);
125:   TSSetEquationType(ts,TS_EQ_IMPLICIT);
126:   TSARKIMEXSetFullyImplicit(ts,PETSC_TRUE);
127:   TSSetIFunction(ts,NULL,(TSIFunction) IFunction,&ctx);
128:   TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,&ctx);

130:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
131:      Set initial conditions
132:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
133:   TSSetSolution(ts,U);

135:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
136:      Set solver options
137:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
138:   TSSetMaxTime(ts,35.0);
139:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);
140:   TSSetTimeStep(ts,.1);
141:   TSSetFromOptions(ts);

143:   direction[0] = direction[1] = 1;
144:   terminate[0] = terminate[1] = PETSC_FALSE;

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

148:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
149:      Solve nonlinear system
150:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
151:   if (ensemble) {
152:     for (du[1] = -2.5; du[1] <= .01; du[1] += .1) {
153:       VecGetArray(U,&u);
154:       u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax);
155:       u[1] = ctx.omega_s;
156:       u[0] += du[0];
157:       u[1] += du[1];
158:       VecRestoreArray(U,&u);
159:       TSSetTimeStep(ts,.01);
160:       TSSolve(ts,U);
161:     }
162:   } else {
163:     TSSolve(ts,U);
164:   }
165:   VecView(U,PETSC_VIEWER_STDOUT_WORLD);
166:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
167:      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
168:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
169:   MatDestroy(&A);
170:   VecDestroy(&U);
171:   TSDestroy(&ts);
172:   PetscFinalize();
173:   return ierr;
174: }


177: /*TEST

179:    build:
180:      requires: !complex !single

182:    test:
183:       args: -nox

185: TEST*/