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