Actual source code: ex1.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{2 H}{\omega_s}\frac{d \omega}{dt} & = & P_m - \frac{EV}{X} \sin(\theta) \\
\frac{d \theta}{dt} = \omega - \omega_s
\end{eqnarray}

 13: /*
 14:    Include "petscts.h" so that we can use TS solvers.  Note that this
 15:    file automatically includes:
 16:      petscsys.h       - base PETSc routines   petscvec.h - vectors
 17:      petscmat.h - matrices
 18:      petscis.h     - index sets            petscksp.h - Krylov subspace methods
 19:      petscviewer.h - viewers               petscpc.h  - preconditioners
 20:      petscksp.h   - linear solvers
 21: */
 22: #include <petscts.h>

 24: typedef struct {
 25:   PetscScalar H,omega_s,E,V,X;
 26:   PetscRandom rand;
 27: } AppCtx;

 29: /*
 30:      Defines the ODE passed to the ODE solver
 31: */
 32: static PetscErrorCode IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,AppCtx *ctx)
 33: {
 34:   PetscErrorCode     ierr;
 35:   PetscScalar        *f,r;
 36:   const PetscScalar  *u,*udot;
 37:   static PetscScalar R = .4;

 40:   PetscRandomGetValue(ctx->rand,&r);
 41:   if (r > .9) R = .5;
 42:   if (r < .1) R = .4;
 43:   R = .4;
 44:   /*  The next three lines allow us to access the entries of the vectors directly */
 45:   VecGetArrayRead(U,&u);
 46:   VecGetArrayRead(Udot,&udot);
 47:   VecGetArray(F,&f);
 48:   f[0] = 2.0*ctx->H*udot[0]/ctx->omega_s + ctx->E*ctx->V*PetscSinScalar(u[1])/ctx->X - R;
 49:   f[1] = udot[1] - u[0] + ctx->omega_s;

 51:   VecRestoreArrayRead(U,&u);
 52:   VecRestoreArrayRead(Udot,&udot);
 53:   VecRestoreArray(F,&f);
 54:   return(0);
 55: }

 57: /*
 58:      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
 59: */
 60: static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,AppCtx *ctx)
 61: {
 62:   PetscErrorCode    ierr;
 63:   PetscInt          rowcol[] = {0,1};
 64:   PetscScalar       J[2][2];
 65:   const PetscScalar *u,*udot;

 68:   VecGetArrayRead(U,&u);
 69:   VecGetArrayRead(Udot,&udot);
 70:   J[0][0] = 2.0*ctx->H*a/ctx->omega_s;   J[0][1] = -ctx->E*ctx->V*PetscCosScalar(u[1])/ctx->X;
 71:   J[1][0] = -1.0;                        J[1][1] = a;
 72:   MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);
 73:   VecRestoreArrayRead(U,&u);
 74:   VecRestoreArrayRead(Udot,&udot);

 76:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 77:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 78:   if (A != B) {
 79:     MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 80:     MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
 81:   }
 82:   return(0);
 83: }

 85: int main(int argc,char **argv)
 86: {
 87:   TS             ts;            /* ODE integrator */
 88:   Vec            U;             /* solution will be stored here */
 89:   Mat            A;             /* Jacobian matrix */
 91:   PetscMPIInt    size;
 92:   PetscInt       n = 2;
 93:   AppCtx         ctx;
 94:   PetscScalar    *u;

 96:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 97:      Initialize program
 98:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 99:   PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
100:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
101:   if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs");

103:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
104:     Create necessary matrix and vectors
105:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
106:   MatCreate(PETSC_COMM_WORLD,&A);
107:   MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);
108:   MatSetFromOptions(A);
109:   MatSetUp(A);

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

113:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
114:     Set runtime options
115:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
116:   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Reaction options","");
117:   {
118:     ctx.omega_s = 1.0;
119:     PetscOptionsScalar("-omega_s","","",ctx.omega_s,&ctx.omega_s,NULL);
120:     ctx.H       = 1.0;
121:     PetscOptionsScalar("-H","","",ctx.H,&ctx.H,NULL);
122:     ctx.E       = 1.0;
123:     PetscOptionsScalar("-E","","",ctx.E,&ctx.E,NULL);
124:     ctx.V       = 1.0;
125:     PetscOptionsScalar("-V","","",ctx.V,&ctx.V,NULL);
126:     ctx.X       = 1.0;
127:     PetscOptionsScalar("-X","","",ctx.X,&ctx.X,NULL);

129:     VecGetArray(U,&u);
130:     u[0] = 1;
131:     u[1] = .7;
132:     VecRestoreArray(U,&u);
133:     PetscOptionsGetVec(NULL,NULL,"-initial",U,NULL);
134:   }
135:   PetscOptionsEnd();

137:   PetscRandomCreate(PETSC_COMM_WORLD,&ctx.rand);
138:   PetscRandomSetFromOptions(ctx.rand);

140:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
141:      Create timestepping solver context
142:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
143:   TSCreate(PETSC_COMM_WORLD,&ts);
144:   TSSetProblemType(ts,TS_NONLINEAR);
145:   TSSetType(ts,TSROSW);
146:   TSSetIFunction(ts,NULL,(TSIFunction) IFunction,&ctx);
147:   TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,&ctx);

149:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
150:      Set initial conditions
151:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
152:   TSSetSolution(ts,U);

154:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
155:      Set solver options
156:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
157:   TSSetMaxTime(ts,2000.0);
158:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);
159:   TSSetTimeStep(ts,.001);
160:   TSSetFromOptions(ts);

162:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
163:      Solve nonlinear system
164:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
165:   TSSolve(ts,U);

167:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
168:      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
169:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
170:   MatDestroy(&A);
171:   VecDestroy(&U);
172:   TSDestroy(&ts);
173:   PetscRandomDestroy(&ctx.rand);
174:   PetscFinalize();
175:   return ierr;
176: }