Actual source code: ex1.c

  1: static char help[] = "An example of hybrid system using TS event.\n";

  3: /*
  4:   The dynamics is described by the ODE
  5:                   u_t = A_i u

  7:   where A_1 = [ 1  -100
  8:                 10  1  ],
  9:         A_2 = [ 1    10
 10:                -100  1 ].
 11:   The index i changes from 1 to 2 when u[1]=2.75u[0] and from 2 to 1 when u[1]=0.36u[0].
 12:   Initially u=[0 1]^T and i=1.

 14:   Reference:
 15:   I. A. Hiskens, M.A. Pai, Trajectory Sensitivity Analysis of Hybrid Systems, IEEE Transactions on Circuits and Systems, Vol 47, No 2, February 2000
 16: */

 18: #include <petscts.h>

 20: typedef struct {
 21:   PetscScalar lambda1;
 22:   PetscScalar lambda2;
 23:   PetscInt    mode; /* mode flag*/
 24: } AppCtx;

 26: PetscErrorCode EventFunction(TS ts, PetscReal t, Vec U, PetscReal *fvalue, void *ctx)
 27: {
 28:   AppCtx            *actx = (AppCtx *)ctx;
 29:   const PetscScalar *u;

 31:   PetscFunctionBegin;
 32:   PetscCall(VecGetArrayRead(U, &u));
 33:   if (actx->mode == 1) {
 34:     fvalue[0] = PetscRealPart(u[1] - actx->lambda1 * u[0]);
 35:   } else if (actx->mode == 2) {
 36:     fvalue[0] = PetscRealPart(u[1] - actx->lambda2 * u[0]);
 37:   }
 38:   PetscCall(VecRestoreArrayRead(U, &u));
 39:   PetscFunctionReturn(PETSC_SUCCESS);
 40: }

 42: PetscErrorCode PostEventFunction(TS ts, PetscInt nevents, PetscInt event_list[], PetscReal t, Vec U, PetscBool forwardsolve, void *ctx)
 43: {
 44:   AppCtx *actx = (AppCtx *)ctx;

 46:   PetscFunctionBegin;
 47:   if (actx->mode == 1) {
 48:     actx->mode = 2;
 49:     PetscCall(PetscPrintf(PETSC_COMM_SELF, "Change from mode 1 to 2 at t = %f \n", (double)t));
 50:   } else if (actx->mode == 2) {
 51:     actx->mode = 1;
 52:     PetscCall(PetscPrintf(PETSC_COMM_SELF, "Change from mode 2 to 1 at t = %f \n", (double)t));
 53:   }
 54:   PetscFunctionReturn(PETSC_SUCCESS);
 55: }

 57: /*
 58:      Defines the ODE passed to the ODE solver
 59: */
 60: static PetscErrorCode IFunction(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, void *ctx)
 61: {
 62:   AppCtx            *actx = (AppCtx *)ctx;
 63:   PetscScalar       *f;
 64:   const PetscScalar *u, *udot;

 66:   PetscFunctionBegin;
 67:   /*  The next three lines allow us to access the entries of the vectors directly */
 68:   PetscCall(VecGetArrayRead(U, &u));
 69:   PetscCall(VecGetArrayRead(Udot, &udot));
 70:   PetscCall(VecGetArray(F, &f));

 72:   if (actx->mode == 1) {
 73:     f[0] = udot[0] - u[0] + 100 * u[1];
 74:     f[1] = udot[1] - 10 * u[0] - u[1];
 75:   } else if (actx->mode == 2) {
 76:     f[0] = udot[0] - u[0] - 10 * u[1];
 77:     f[1] = udot[1] + 100 * u[0] - u[1];
 78:   }

 80:   PetscCall(VecRestoreArrayRead(U, &u));
 81:   PetscCall(VecRestoreArrayRead(Udot, &udot));
 82:   PetscCall(VecRestoreArray(F, &f));
 83:   PetscFunctionReturn(PETSC_SUCCESS);
 84: }

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

 96:   PetscFunctionBegin;
 97:   PetscCall(VecGetArrayRead(U, &u));
 98:   PetscCall(VecGetArrayRead(Udot, &udot));

100:   if (actx->mode == 1) {
101:     J[0][0] = a - 1;
102:     J[0][1] = 100;
103:     J[1][0] = -10;
104:     J[1][1] = a - 1;
105:   } else if (actx->mode == 2) {
106:     J[0][0] = a - 1;
107:     J[0][1] = -10;
108:     J[1][0] = 100;
109:     J[1][1] = a - 1;
110:   }
111:   PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES));

113:   PetscCall(VecRestoreArrayRead(U, &u));
114:   PetscCall(VecRestoreArrayRead(Udot, &udot));

116:   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
117:   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
118:   if (A != B) {
119:     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
120:     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
121:   }
122:   PetscFunctionReturn(PETSC_SUCCESS);
123: }

125: int main(int argc, char **argv)
126: {
127:   TS           ts; /* ODE integrator */
128:   Vec          U;  /* solution will be stored here */
129:   Mat          A;  /* Jacobian matrix */
130:   PetscMPIInt  size;
131:   PetscInt     n = 2;
132:   PetscScalar *u;
133:   AppCtx       app;
134:   PetscInt     direction[1];
135:   PetscBool    terminate[1];

137:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
138:      Initialize program
139:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
140:   PetscFunctionBeginUser;
141:   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
142:   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
143:   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs");
144:   app.mode    = 1;
145:   app.lambda1 = 2.75;
146:   app.lambda2 = 0.36;
147:   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "ex1 options", "");
148:   {
149:     PetscCall(PetscOptionsReal("-lambda1", "", "", app.lambda1, &app.lambda1, NULL));
150:     PetscCall(PetscOptionsReal("-lambda2", "", "", app.lambda2, &app.lambda2, NULL));
151:   }
152:   PetscOptionsEnd();

154:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
155:     Create necessary matrix and vectors
156:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
157:   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
158:   PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE));
159:   PetscCall(MatSetType(A, MATDENSE));
160:   PetscCall(MatSetFromOptions(A));
161:   PetscCall(MatSetUp(A));

163:   PetscCall(MatCreateVecs(A, &U, NULL));

165:   PetscCall(VecGetArray(U, &u));
166:   u[0] = 0;
167:   u[1] = 1;
168:   PetscCall(VecRestoreArray(U, &u));
169:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
170:      Create timestepping solver context
171:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
172:   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
173:   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
174:   PetscCall(TSSetType(ts, TSCN));
175:   PetscCall(TSSetIFunction(ts, NULL, (TSIFunctionFn *)IFunction, &app));
176:   PetscCall(TSSetIJacobian(ts, A, A, (TSIJacobianFn *)IJacobian, &app));

178:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
179:      Set initial conditions
180:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
181:   PetscCall(TSSetSolution(ts, U));

183:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
184:      Set solver options
185:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
186:   PetscCall(TSSetMaxTime(ts, 0.125));
187:   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
188:   PetscCall(TSSetTimeStep(ts, 1. / 256.));
189:   PetscCall(TSSetFromOptions(ts));

191:   /* Set directions and terminate flags for the two events */
192:   direction[0] = 0;
193:   terminate[0] = PETSC_FALSE;
194:   PetscCall(TSSetEventHandler(ts, 1, direction, terminate, EventFunction, PostEventFunction, (void *)&app));

196:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
197:      Run timestepping solver
198:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
199:   PetscCall(TSSolve(ts, U));

201:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
202:      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
203:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
204:   PetscCall(MatDestroy(&A));
205:   PetscCall(VecDestroy(&U));
206:   PetscCall(TSDestroy(&ts));

208:   PetscCall(PetscFinalize());
209:   return 0;
210: }

212: /*TEST

214:    build:
215:       requires: !complex
216:    test:
217:       args: -ts_monitor

219:    test:
220:       suffix: 2
221:       args: -ts_monitor_lg_solution -1
222:       requires: x

224: TEST*/