1 static char help[] = "An example of hybrid system using TS event.\n"; 2 3 /* 4 The dynamics is described by the ODE 5 u_t = A_i u 6 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. 13 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 */ 17 18 #include <petscts.h> 19 20 typedef struct { 21 PetscScalar lambda1; 22 PetscScalar lambda2; 23 PetscInt mode; /* mode flag*/ 24 } AppCtx; 25 26 PetscErrorCode EventFunction(TS ts,PetscReal t,Vec U,PetscScalar *fvalue,void *ctx) 27 { 28 AppCtx *actx=(AppCtx*)ctx; 29 PetscErrorCode ierr; 30 const PetscScalar *u; 31 32 PetscFunctionBegin; 33 ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 34 if (actx->mode == 1) { 35 fvalue[0] = u[1]-actx->lambda1*u[0]; 36 }else if (actx->mode == 2) { 37 fvalue[0] = u[1]-actx->lambda2*u[0]; 38 } 39 ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 40 PetscFunctionReturn(0); 41 } 42 43 PetscErrorCode PostEventFunction(TS ts,PetscInt nevents,PetscInt event_list[],PetscReal t,Vec U,PetscBool forwardsolve,void* ctx) 44 { 45 AppCtx *actx=(AppCtx*)ctx; 46 PetscErrorCode ierr; 47 48 PetscFunctionBegin; 49 if (actx->mode == 1) { 50 actx->mode = 2; 51 ierr = PetscPrintf(PETSC_COMM_SELF,"Change from mode 1 to 2 at t = %f \n",(double)t);CHKERRQ(ierr); 52 } else if (actx->mode == 2) { 53 actx->mode = 1; 54 ierr = PetscPrintf(PETSC_COMM_SELF,"Change from mode 2 to 1 at t = %f \n",(double)t);CHKERRQ(ierr); 55 } 56 PetscFunctionReturn(0); 57 } 58 59 /* 60 Defines the ODE passed to the ODE solver 61 */ 62 static PetscErrorCode IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,void *ctx) 63 { 64 AppCtx *actx=(AppCtx*)ctx; 65 PetscErrorCode ierr; 66 PetscScalar *f; 67 const PetscScalar *u,*udot; 68 69 PetscFunctionBegin; 70 /* The next three lines allow us to access the entries of the vectors directly */ 71 ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 72 ierr = VecGetArrayRead(Udot,&udot);CHKERRQ(ierr); 73 ierr = VecGetArray(F,&f);CHKERRQ(ierr); 74 75 if (actx->mode == 1) { 76 f[0] = udot[0]-u[0]+100*u[1]; 77 f[1] = udot[1]-10*u[0]-u[1]; 78 } else if (actx->mode == 2) { 79 f[0] = udot[0]-u[0]-10*u[1]; 80 f[1] = udot[1]+100*u[0]-u[1]; 81 } 82 83 ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 84 ierr = VecRestoreArrayRead(Udot,&udot);CHKERRQ(ierr); 85 ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 86 PetscFunctionReturn(0); 87 } 88 89 /* 90 Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian. 91 */ 92 static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,void *ctx) 93 { 94 AppCtx *actx=(AppCtx*)ctx; 95 PetscErrorCode ierr; 96 PetscInt rowcol[] = {0,1}; 97 PetscScalar J[2][2]; 98 const PetscScalar *u,*udot; 99 100 PetscFunctionBegin; 101 ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 102 ierr = VecGetArrayRead(Udot,&udot);CHKERRQ(ierr); 103 104 if (actx->mode == 1) { 105 J[0][0] = a-1; J[0][1] = 100; 106 J[1][0] = -10; J[1][1] = a-1; 107 } else if (actx->mode == 2) { 108 J[0][0] = a-1; J[0][1] = -10; 109 J[1][0] = 100; J[1][1] = a-1; 110 } 111 ierr = MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);CHKERRQ(ierr); 112 113 ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 114 ierr = VecRestoreArrayRead(Udot,&udot);CHKERRQ(ierr); 115 116 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 117 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 118 if (A != B) { 119 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 120 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 121 } 122 PetscFunctionReturn(0); 123 } 124 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 PetscErrorCode ierr; 131 PetscMPIInt size; 132 PetscInt n = 2; 133 PetscScalar *u; 134 AppCtx app; 135 PetscInt direction[1]; 136 PetscBool terminate[1]; 137 138 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 139 Initialize program 140 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 141 ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 142 ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 143 if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs"); 144 app.mode = 1; 145 app.lambda1 = 2.75; 146 app.lambda2 = 0.36; 147 ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"ex1 options","");CHKERRQ(ierr); 148 { 149 ierr = PetscOptionsReal("-lambda1","","",app.lambda1,&app.lambda1,NULL);CHKERRQ(ierr); 150 ierr = PetscOptionsReal("-lambda2","","",app.lambda2,&app.lambda2,NULL);CHKERRQ(ierr); 151 } 152 ierr = PetscOptionsEnd();CHKERRQ(ierr); 153 154 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 155 Create necessary matrix and vectors 156 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 157 ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 158 ierr = MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 159 ierr = MatSetType(A,MATDENSE);CHKERRQ(ierr); 160 ierr = MatSetFromOptions(A);CHKERRQ(ierr); 161 ierr = MatSetUp(A);CHKERRQ(ierr); 162 163 ierr = MatCreateVecs(A,&U,NULL);CHKERRQ(ierr); 164 165 ierr = VecGetArray(U,&u);CHKERRQ(ierr); 166 u[0] = 0; 167 u[1] = 1; 168 ierr = VecRestoreArray(U,&u);CHKERRQ(ierr); 169 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 170 Create timestepping solver context 171 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 172 ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr); 173 ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr); 174 ierr = TSSetType(ts,TSCN);CHKERRQ(ierr); 175 ierr = TSSetIFunction(ts,NULL,(TSIFunction)IFunction,&app);CHKERRQ(ierr); 176 ierr = TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,&app);CHKERRQ(ierr); 177 178 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 179 Set initial conditions 180 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 181 ierr = TSSetSolution(ts,U);CHKERRQ(ierr); 182 183 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 184 Set solver options 185 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 186 ierr = TSSetMaxTime(ts,0.125);CHKERRQ(ierr); 187 ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr); 188 ierr = TSSetTimeStep(ts,1./256.);CHKERRQ(ierr); 189 ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 190 191 /* Set directions and terminate flags for the two events */ 192 direction[0] = 0; 193 terminate[0] = PETSC_FALSE; 194 ierr = TSSetEventHandler(ts,1,direction,terminate,EventFunction,PostEventFunction,(void*)&app);CHKERRQ(ierr); 195 196 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 197 Run timestepping solver 198 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 199 ierr = TSSolve(ts,U);CHKERRQ(ierr); 200 201 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 202 Free work space. All PETSc objects should be destroyed when they are no longer needed. 203 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 204 ierr = MatDestroy(&A);CHKERRQ(ierr); 205 ierr = VecDestroy(&U);CHKERRQ(ierr); 206 ierr = TSDestroy(&ts);CHKERRQ(ierr); 207 208 ierr = PetscFinalize(); 209 return(ierr); 210 } 211 212 213 /*TEST 214 215 build: 216 requires: !complex 217 test: 218 args: -ts_monitor 219 220 test: 221 suffix: 2 222 args: -ts_monitor_lg_solution -1 223 requires: x 224 225 TEST*/ 226