1c4762a1bSJed Brown static char help[] = "Trajectory sensitivity of a hybrid system with state-dependent switchings.\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 References: 15*606c0280SSatish Balay + * - H. Zhang, S. Abhyankar, E. Constantinescu, M. Mihai, Discrete Adjoint Sensitivity Analysis of Hybrid Dynamical Systems With Switching, IEEE Transactions on Circuits and Systems I: Regular Papers, 64(5), May 2017 16*606c0280SSatish Balay - * - I. A. Hiskens, M.A. Pai, Trajectory Sensitivity Analysis of Hybrid Systems, IEEE Transactions on Circuits and Systems, Vol 47, No 2, February 2000 17c4762a1bSJed Brown */ 18c4762a1bSJed Brown 19c4762a1bSJed Brown #include <petscts.h> 20c4762a1bSJed Brown 21c4762a1bSJed Brown typedef struct { 22c4762a1bSJed Brown PetscScalar lambda1; 23c4762a1bSJed Brown PetscScalar lambda2; 24c4762a1bSJed Brown PetscInt mode; /* mode flag*/ 25c4762a1bSJed Brown PetscReal print_time; 26c4762a1bSJed Brown } AppCtx; 27c4762a1bSJed Brown 28c4762a1bSJed Brown PetscErrorCode MyMonitor(TS ts,PetscInt stepnum,PetscReal time,Vec U,void *ctx) 29c4762a1bSJed Brown { 30c4762a1bSJed Brown AppCtx *actx=(AppCtx*)ctx; 31c4762a1bSJed Brown Mat sp; 32c4762a1bSJed Brown PetscScalar *u; 33c4762a1bSJed Brown PetscInt nump; 34c4762a1bSJed Brown FILE *f; 35c4762a1bSJed Brown PetscErrorCode ierr; 36c4762a1bSJed Brown 37c4762a1bSJed Brown PetscFunctionBegin; 38c4762a1bSJed Brown if (time >= actx->print_time) { 39c4762a1bSJed Brown actx->print_time += 1./256.; 40c4762a1bSJed Brown ierr = TSForwardGetSensitivities(ts,&nump,&sp);CHKERRQ(ierr); 41c4762a1bSJed Brown ierr = MatDenseGetColumn(sp,2,&u);CHKERRQ(ierr); 42c4762a1bSJed Brown f = fopen("fwd_sp.out", "a"); 43c4762a1bSJed Brown ierr = PetscFPrintf(PETSC_COMM_WORLD,f,"%20.15lf %20.15lf %20.15lf\n",time,u[0],u[1]);CHKERRQ(ierr); 44c4762a1bSJed Brown ierr = MatDenseRestoreColumn(sp,&u);CHKERRQ(ierr); 45c4762a1bSJed Brown fclose(f); 46c4762a1bSJed Brown } 47c4762a1bSJed Brown PetscFunctionReturn(0); 48c4762a1bSJed Brown } 49c4762a1bSJed Brown 50c4762a1bSJed Brown PetscErrorCode EventFunction(TS ts,PetscReal t,Vec U,PetscScalar *fvalue,void *ctx) 51c4762a1bSJed Brown { 52c4762a1bSJed Brown AppCtx *actx=(AppCtx*)ctx; 53c4762a1bSJed Brown PetscErrorCode ierr; 54c4762a1bSJed Brown const PetscScalar *u; 55c4762a1bSJed Brown 56c4762a1bSJed Brown PetscFunctionBegin; 57c4762a1bSJed Brown ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 58c4762a1bSJed Brown if (actx->mode == 1) { 59c4762a1bSJed Brown fvalue[0] = u[1]-actx->lambda1*u[0]; 60c4762a1bSJed Brown }else if (actx->mode == 2) { 61c4762a1bSJed Brown fvalue[0] = u[1]-actx->lambda2*u[0]; 62c4762a1bSJed Brown } 63c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 64c4762a1bSJed Brown PetscFunctionReturn(0); 65c4762a1bSJed Brown } 66c4762a1bSJed Brown 67c4762a1bSJed Brown PetscErrorCode ShiftGradients(TS ts,Vec U,AppCtx *actx) 68c4762a1bSJed Brown { 69c4762a1bSJed Brown Mat sp; 70c4762a1bSJed Brown PetscScalar *x; 71c4762a1bSJed Brown PetscScalar *u; 72c4762a1bSJed Brown PetscErrorCode ierr; 73c4762a1bSJed Brown PetscScalar tmp[2],A1[2][2],A2[2],denorm; 74c4762a1bSJed Brown PetscInt nump; 75c4762a1bSJed Brown 76c4762a1bSJed Brown PetscFunctionBegin; 77c4762a1bSJed Brown ierr = TSForwardGetSensitivities(ts,&nump,&sp);CHKERRQ(ierr); 78c4762a1bSJed Brown ierr = VecGetArray(U,&u);CHKERRQ(ierr); 79c4762a1bSJed Brown 80c4762a1bSJed Brown if (actx->mode==1) { 81c4762a1bSJed Brown denorm = -actx->lambda1*(u[0]-100.*u[1])+1.*(10.*u[0]+u[1]); 82c4762a1bSJed Brown A1[0][0] = 110.*u[1]*(-actx->lambda1)/denorm+1.; 83c4762a1bSJed Brown A1[1][0] = -110.*u[0]*(-actx->lambda1)/denorm; 84c4762a1bSJed Brown A1[0][1] = 110.*u[1]*1./denorm; 85c4762a1bSJed Brown A1[1][1] = -110.*u[0]*1./denorm+1.; 86c4762a1bSJed Brown 87c4762a1bSJed Brown A2[0] = 110.*u[1]*(-u[0])/denorm; 88c4762a1bSJed Brown A2[1] = -110.*u[0]*(-u[0])/denorm; 89c4762a1bSJed Brown }else { 90c4762a1bSJed Brown denorm = -actx->lambda2*(u[0]+10.*u[1])+1.*(-100.*u[0]+u[1]); 91c4762a1bSJed Brown A1[0][0] = 110.*u[1]*(actx->lambda2)/denorm+1; 92c4762a1bSJed Brown A1[1][0] = -110.*u[0]*(actx->lambda2)/denorm; 93c4762a1bSJed Brown A1[0][1] = -110.*u[1]*1./denorm; 94c4762a1bSJed Brown A1[1][1] = 110.*u[0]*1./denorm+1.; 95c4762a1bSJed Brown 96c4762a1bSJed Brown A2[0] = 0; 97c4762a1bSJed Brown A2[1] = 0; 98c4762a1bSJed Brown } 99c4762a1bSJed Brown 100c4762a1bSJed Brown ierr = VecRestoreArray(U,&u);CHKERRQ(ierr); 101c4762a1bSJed Brown 102c4762a1bSJed Brown ierr = MatDenseGetColumn(sp,0,&x);CHKERRQ(ierr); 103c4762a1bSJed Brown tmp[0] = A1[0][0]*x[0]+A1[0][1]*x[1]; 104c4762a1bSJed Brown tmp[1] = A1[1][0]*x[0]+A1[1][1]*x[1]; 105c4762a1bSJed Brown x[0] = tmp[0]; 106c4762a1bSJed Brown x[1] = tmp[1]; 107c4762a1bSJed Brown ierr = MatDenseRestoreColumn(sp,&x);CHKERRQ(ierr); 108c4762a1bSJed Brown 109c4762a1bSJed Brown ierr = MatDenseGetColumn(sp,1,&x);CHKERRQ(ierr); 110c4762a1bSJed Brown tmp[0] = A1[0][0]*x[0]+A1[0][1]*x[1]; 111c4762a1bSJed Brown tmp[1] = A1[1][0]*x[0]+A1[1][1]*x[1]; 112c4762a1bSJed Brown x[0] = tmp[0]; 113c4762a1bSJed Brown x[1] = tmp[1]; 114c4762a1bSJed Brown ierr = MatDenseRestoreColumn(sp,&x);CHKERRQ(ierr); 115c4762a1bSJed Brown 116c4762a1bSJed Brown ierr = MatDenseGetColumn(sp,2,&x);CHKERRQ(ierr); 117c4762a1bSJed Brown tmp[0] = A1[0][0]*x[0]+A1[0][1]*x[1]; 118c4762a1bSJed Brown tmp[1] = A1[1][0]*x[0]+A1[1][1]*x[1]; 119c4762a1bSJed Brown x[0] = tmp[0]+A2[0]; 120c4762a1bSJed Brown x[1] = tmp[1]+A2[1]; 121c4762a1bSJed Brown ierr = MatDenseRestoreColumn(sp,&x);CHKERRQ(ierr); 122c4762a1bSJed Brown 123c4762a1bSJed Brown PetscFunctionReturn(0); 124c4762a1bSJed Brown } 125c4762a1bSJed Brown 126c4762a1bSJed Brown PetscErrorCode PostEventFunction(TS ts,PetscInt nevents,PetscInt event_list[],PetscReal t,Vec U,PetscBool forwardsolve,void* ctx) 127c4762a1bSJed Brown { 128c4762a1bSJed Brown AppCtx *actx=(AppCtx*)ctx; 129c4762a1bSJed Brown PetscErrorCode ierr; 130c4762a1bSJed Brown 131c4762a1bSJed Brown PetscFunctionBegin; 132c4762a1bSJed Brown /* ierr = VecView(U,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */ 133c4762a1bSJed Brown ierr = ShiftGradients(ts,U,actx);CHKERRQ(ierr); 134c4762a1bSJed Brown if (actx->mode == 1) { 135c4762a1bSJed Brown actx->mode = 2; 136c4762a1bSJed Brown /* ierr = PetscPrintf(PETSC_COMM_SELF,"Change from mode 1 to 2 at t = %f \n",(double)t);CHKERRQ(ierr); */ 137c4762a1bSJed Brown } else if (actx->mode == 2) { 138c4762a1bSJed Brown actx->mode = 1; 139c4762a1bSJed Brown /* ierr = PetscPrintf(PETSC_COMM_SELF,"Change from mode 2 to 1 at t = %f \n",(double)t);CHKERRQ(ierr); */ 140c4762a1bSJed Brown } 141c4762a1bSJed Brown PetscFunctionReturn(0); 142c4762a1bSJed Brown } 143c4762a1bSJed Brown 144c4762a1bSJed Brown /* 145c4762a1bSJed Brown Defines the ODE passed to the ODE solver 146c4762a1bSJed Brown */ 147c4762a1bSJed Brown static PetscErrorCode IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,void *ctx) 148c4762a1bSJed Brown { 149c4762a1bSJed Brown AppCtx *actx=(AppCtx*)ctx; 150c4762a1bSJed Brown PetscErrorCode ierr; 151c4762a1bSJed Brown PetscScalar *f; 152c4762a1bSJed Brown const PetscScalar *u,*udot; 153c4762a1bSJed Brown 154c4762a1bSJed Brown PetscFunctionBegin; 155c4762a1bSJed Brown /* The next three lines allow us to access the entries of the vectors directly */ 156c4762a1bSJed Brown ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 157c4762a1bSJed Brown ierr = VecGetArrayRead(Udot,&udot);CHKERRQ(ierr); 158c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 159c4762a1bSJed Brown 160c4762a1bSJed Brown if (actx->mode == 1) { 161c4762a1bSJed Brown f[0] = udot[0]-u[0]+100*u[1]; 162c4762a1bSJed Brown f[1] = udot[1]-10*u[0]-u[1]; 163c4762a1bSJed Brown } else if (actx->mode == 2) { 164c4762a1bSJed Brown f[0] = udot[0]-u[0]-10*u[1]; 165c4762a1bSJed Brown f[1] = udot[1]+100*u[0]-u[1]; 166c4762a1bSJed Brown } 167c4762a1bSJed Brown 168c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 169c4762a1bSJed Brown ierr = VecRestoreArrayRead(Udot,&udot);CHKERRQ(ierr); 170c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 171c4762a1bSJed Brown PetscFunctionReturn(0); 172c4762a1bSJed Brown } 173c4762a1bSJed Brown 174c4762a1bSJed Brown /* 175c4762a1bSJed Brown Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian. 176c4762a1bSJed Brown */ 177c4762a1bSJed Brown static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,void *ctx) 178c4762a1bSJed Brown { 179c4762a1bSJed Brown AppCtx *actx=(AppCtx*)ctx; 180c4762a1bSJed Brown PetscErrorCode ierr; 181c4762a1bSJed Brown PetscInt rowcol[] = {0,1}; 182c4762a1bSJed Brown PetscScalar J[2][2]; 183c4762a1bSJed Brown const PetscScalar *u,*udot; 184c4762a1bSJed Brown 185c4762a1bSJed Brown PetscFunctionBegin; 186c4762a1bSJed Brown ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 187c4762a1bSJed Brown ierr = VecGetArrayRead(Udot,&udot);CHKERRQ(ierr); 188c4762a1bSJed Brown 189c4762a1bSJed Brown if (actx->mode == 1) { 190c4762a1bSJed Brown J[0][0] = a-1; J[0][1] = 100; 191c4762a1bSJed Brown J[1][0] = -10; J[1][1] = a-1; 192c4762a1bSJed Brown } else if (actx->mode == 2) { 193c4762a1bSJed Brown J[0][0] = a-1; J[0][1] = -10; 194c4762a1bSJed Brown J[1][0] = 100; J[1][1] = a-1; 195c4762a1bSJed Brown } 196c4762a1bSJed Brown ierr = MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);CHKERRQ(ierr); 197c4762a1bSJed Brown 198c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 199c4762a1bSJed Brown ierr = VecRestoreArrayRead(Udot,&udot);CHKERRQ(ierr); 200c4762a1bSJed Brown 201c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 202c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 203c4762a1bSJed Brown if (A != B) { 204c4762a1bSJed Brown ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 205c4762a1bSJed Brown ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 206c4762a1bSJed Brown } 207c4762a1bSJed Brown PetscFunctionReturn(0); 208c4762a1bSJed Brown } 209c4762a1bSJed Brown 210c4762a1bSJed Brown /* Matrix JacobianP is constant so that it only needs to be evaluated once */ 211c4762a1bSJed Brown static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat Ap,void *ctx) 212c4762a1bSJed Brown { 213c4762a1bSJed Brown PetscFunctionBeginUser; 214c4762a1bSJed Brown PetscFunctionReturn(0); 215c4762a1bSJed Brown } 216c4762a1bSJed Brown 217c4762a1bSJed Brown int main(int argc,char **argv) 218c4762a1bSJed Brown { 219c4762a1bSJed Brown TS ts; /* ODE integrator */ 220c4762a1bSJed Brown Vec U; /* solution will be stored here */ 221c4762a1bSJed Brown Mat A; /* Jacobian matrix */ 222c4762a1bSJed Brown Mat Ap; /* Jacobian dfdp */ 223c4762a1bSJed Brown PetscErrorCode ierr; 224c4762a1bSJed Brown PetscMPIInt size; 225c4762a1bSJed Brown PetscInt n = 2; 226c4762a1bSJed Brown PetscScalar *u; 227c4762a1bSJed Brown AppCtx app; 228c4762a1bSJed Brown PetscInt direction[1]; 229c4762a1bSJed Brown PetscBool terminate[1]; 230c4762a1bSJed Brown Mat sp; 231c4762a1bSJed Brown PetscReal tend; 232c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 233c4762a1bSJed Brown Initialize program 234c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 235c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 236ffc4695bSBarry Smith ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 2372c71b3e2SJacob Faibussowitsch PetscCheckFalse(size > 1,PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs"); 238c4762a1bSJed Brown app.mode = 1; 239c4762a1bSJed Brown app.lambda1 = 2.75; 240c4762a1bSJed Brown app.lambda2 = 0.36; 241c4762a1bSJed Brown app.print_time = 1./256.; 242c4762a1bSJed Brown tend = 0.125; 243c4762a1bSJed Brown ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"ex1fwd options","");CHKERRQ(ierr); 244c4762a1bSJed Brown { 245c4762a1bSJed Brown ierr = PetscOptionsReal("-lambda1","","",app.lambda1,&app.lambda1,NULL);CHKERRQ(ierr); 246c4762a1bSJed Brown ierr = PetscOptionsReal("-lambda2","","",app.lambda2,&app.lambda2,NULL);CHKERRQ(ierr); 247c4762a1bSJed Brown ierr = PetscOptionsReal("-tend","","",tend,&tend,NULL);CHKERRQ(ierr); 248c4762a1bSJed Brown } 249c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 250c4762a1bSJed Brown 251c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 252c4762a1bSJed Brown Create necessary matrix and vectors 253c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 254c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 255c4762a1bSJed Brown ierr = MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 256c4762a1bSJed Brown ierr = MatSetType(A,MATDENSE);CHKERRQ(ierr); 257c4762a1bSJed Brown ierr = MatSetFromOptions(A);CHKERRQ(ierr); 258c4762a1bSJed Brown ierr = MatSetUp(A);CHKERRQ(ierr); 259c4762a1bSJed Brown 260c4762a1bSJed Brown ierr = MatCreateVecs(A,&U,NULL);CHKERRQ(ierr); 261c4762a1bSJed Brown 262c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&Ap);CHKERRQ(ierr); 263c4762a1bSJed Brown ierr = MatSetSizes(Ap,n,3,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 264c4762a1bSJed Brown ierr = MatSetType(Ap,MATDENSE);CHKERRQ(ierr); 265c4762a1bSJed Brown ierr = MatSetFromOptions(Ap);CHKERRQ(ierr); 266c4762a1bSJed Brown ierr = MatSetUp(Ap);CHKERRQ(ierr); 267c4762a1bSJed Brown ierr = MatZeroEntries(Ap);CHKERRQ(ierr); 268c4762a1bSJed Brown 269c4762a1bSJed Brown ierr = MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,3,NULL,&sp);CHKERRQ(ierr); 270c4762a1bSJed Brown ierr = MatZeroEntries(sp);CHKERRQ(ierr); 271c4762a1bSJed Brown ierr = MatShift(sp,1.0);CHKERRQ(ierr); 272c4762a1bSJed Brown 273c4762a1bSJed Brown ierr = VecGetArray(U,&u);CHKERRQ(ierr); 274c4762a1bSJed Brown u[0] = 0; 275c4762a1bSJed Brown u[1] = 1; 276c4762a1bSJed Brown ierr = VecRestoreArray(U,&u);CHKERRQ(ierr); 277c4762a1bSJed Brown 278c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 279c4762a1bSJed Brown Create timestepping solver context 280c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 281c4762a1bSJed Brown ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr); 282c4762a1bSJed Brown ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr); 283c4762a1bSJed Brown ierr = TSSetType(ts,TSCN);CHKERRQ(ierr); 284c4762a1bSJed Brown ierr = TSSetIFunction(ts,NULL,(TSIFunction)IFunction,&app);CHKERRQ(ierr); 285c4762a1bSJed Brown ierr = TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,&app);CHKERRQ(ierr); 286c4762a1bSJed Brown 287c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 288c4762a1bSJed Brown Set initial conditions 289c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 290c4762a1bSJed Brown ierr = TSSetSolution(ts,U);CHKERRQ(ierr); 291c4762a1bSJed Brown 292c4762a1bSJed Brown ierr = TSForwardSetSensitivities(ts,3,sp);CHKERRQ(ierr); 293c4762a1bSJed Brown /* Set RHS JacobianP */ 294c4762a1bSJed Brown ierr = TSSetRHSJacobianP(ts,Ap,RHSJacobianP,&app);CHKERRQ(ierr); 295c4762a1bSJed Brown 296c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 297c4762a1bSJed Brown Set solver options 298c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 299c4762a1bSJed Brown ierr = TSSetMaxTime(ts,tend);CHKERRQ(ierr); 300c4762a1bSJed Brown ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr); 301c4762a1bSJed Brown ierr = TSSetTimeStep(ts,1./256.);CHKERRQ(ierr); 302c4762a1bSJed Brown ierr = TSMonitorSet(ts,MyMonitor,&app,PETSC_NULL);CHKERRQ(ierr); 303c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 304c4762a1bSJed Brown 305c4762a1bSJed Brown /* Set directions and terminate flags for the two events */ 306c4762a1bSJed Brown direction[0] = 0; 307c4762a1bSJed Brown terminate[0] = PETSC_FALSE; 308c4762a1bSJed Brown ierr = TSSetEventHandler(ts,1,direction,terminate,EventFunction,PostEventFunction,(void*)&app);CHKERRQ(ierr); 309c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 310c4762a1bSJed Brown Run timestepping solver 311c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 312c4762a1bSJed Brown ierr = TSSolve(ts,U);CHKERRQ(ierr); 313c4762a1bSJed Brown 314c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 315c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they are no longer needed. 316c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 317c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 318c4762a1bSJed Brown ierr = VecDestroy(&U);CHKERRQ(ierr); 319c4762a1bSJed Brown ierr = TSDestroy(&ts);CHKERRQ(ierr); 320c4762a1bSJed Brown 321c4762a1bSJed Brown ierr = MatDestroy(&Ap);CHKERRQ(ierr); 322c4762a1bSJed Brown ierr = MatDestroy(&sp);CHKERRQ(ierr); 323c4762a1bSJed Brown ierr = PetscFinalize(); 324c4762a1bSJed Brown return(ierr); 325c4762a1bSJed Brown } 326c4762a1bSJed Brown 327c4762a1bSJed Brown /*TEST 328c4762a1bSJed Brown 329c4762a1bSJed Brown build: 330c4762a1bSJed Brown requires: !complex 331c4762a1bSJed Brown 332c4762a1bSJed Brown test: 333c4762a1bSJed Brown args: -ts_monitor 334c4762a1bSJed Brown 335c4762a1bSJed Brown TEST*/ 336