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