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 PetscReal lambda1; 22*c4762a1bSJed Brown PetscReal lambda2; 23*c4762a1bSJed Brown PetscInt mode; /* mode flag*/ 24*c4762a1bSJed Brown } AppCtx; 25*c4762a1bSJed Brown 26*c4762a1bSJed Brown PetscErrorCode FWDRun(TS, Vec, void *); 27*c4762a1bSJed Brown 28*c4762a1bSJed Brown PetscErrorCode EventFunction(TS ts,PetscReal t,Vec U,PetscScalar *fvalue,void *ctx) 29*c4762a1bSJed Brown { 30*c4762a1bSJed Brown AppCtx *actx=(AppCtx*)ctx; 31*c4762a1bSJed Brown PetscErrorCode ierr; 32*c4762a1bSJed Brown const PetscScalar *u; 33*c4762a1bSJed Brown 34*c4762a1bSJed Brown PetscFunctionBegin; 35*c4762a1bSJed Brown ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 36*c4762a1bSJed Brown if (actx->mode == 1) { 37*c4762a1bSJed Brown fvalue[0] = u[1]-actx->lambda1*u[0]; 38*c4762a1bSJed Brown } else if (actx->mode == 2) { 39*c4762a1bSJed Brown fvalue[0] = u[1]-actx->lambda2*u[0]; 40*c4762a1bSJed Brown } 41*c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 42*c4762a1bSJed Brown PetscFunctionReturn(0); 43*c4762a1bSJed Brown } 44*c4762a1bSJed Brown 45*c4762a1bSJed Brown PetscErrorCode ShiftGradients(TS ts,Vec U,AppCtx *actx) 46*c4762a1bSJed Brown { 47*c4762a1bSJed Brown Vec *lambda,*mu; 48*c4762a1bSJed Brown PetscScalar *x,*y; 49*c4762a1bSJed Brown const PetscScalar *u; 50*c4762a1bSJed Brown PetscErrorCode ierr; 51*c4762a1bSJed Brown PetscScalar tmp[2],A1[2][2],A2[2],denorm1,denorm2; 52*c4762a1bSJed Brown PetscInt numcost; 53*c4762a1bSJed Brown 54*c4762a1bSJed Brown PetscFunctionBegin; 55*c4762a1bSJed Brown ierr = TSGetCostGradients(ts,&numcost,&lambda,&mu);CHKERRQ(ierr); 56*c4762a1bSJed Brown ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 57*c4762a1bSJed Brown 58*c4762a1bSJed Brown if (actx->mode==2) { 59*c4762a1bSJed Brown denorm1 = -actx->lambda1*(u[0]-100.*u[1])+1.*(10.*u[0]+u[1]); 60*c4762a1bSJed Brown denorm2 = -actx->lambda1*(u[0]+10.*u[1])+1.*(-100.*u[0]+u[1]); 61*c4762a1bSJed Brown A1[0][0] = 110.*u[1]*(-actx->lambda1)/denorm1+1.; 62*c4762a1bSJed Brown A1[0][1] = -110.*u[0]*(-actx->lambda1)/denorm1; 63*c4762a1bSJed Brown A1[1][0] = 110.*u[1]*1./denorm1; 64*c4762a1bSJed Brown A1[1][1] = -110.*u[0]*1./denorm1+1.; 65*c4762a1bSJed Brown 66*c4762a1bSJed Brown A2[0] = 110.*u[1]*(-u[0])/denorm2; 67*c4762a1bSJed Brown A2[1] = -110.*u[0]*(-u[0])/denorm2; 68*c4762a1bSJed Brown } else { 69*c4762a1bSJed Brown denorm2 = -actx->lambda2*(u[0]+10.*u[1])+1.*(-100.*u[0]+u[1]); 70*c4762a1bSJed Brown A1[0][0] = 110.*u[1]*(-actx->lambda1)/denorm2+1; 71*c4762a1bSJed Brown A1[0][1] = -110.*u[0]*(-actx->lambda1)/denorm2; 72*c4762a1bSJed Brown A1[1][0] = 110.*u[1]*1./denorm2; 73*c4762a1bSJed Brown A1[1][1] = -110.*u[0]*1./denorm2+1.; 74*c4762a1bSJed Brown 75*c4762a1bSJed Brown A2[0] = 0; 76*c4762a1bSJed Brown A2[1] = 0; 77*c4762a1bSJed Brown } 78*c4762a1bSJed Brown 79*c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 80*c4762a1bSJed Brown 81*c4762a1bSJed Brown ierr = VecGetArray(lambda[0],&x);CHKERRQ(ierr); 82*c4762a1bSJed Brown ierr = VecGetArray(mu[0],&y);CHKERRQ(ierr); 83*c4762a1bSJed Brown tmp[0] = A1[0][0]*x[0]+A1[0][1]*x[1]; 84*c4762a1bSJed Brown tmp[1] = A1[1][0]*x[0]+A1[1][1]*x[1]; 85*c4762a1bSJed Brown y[0] = y[0] + A2[0]*x[0]+A2[1]*x[1]; 86*c4762a1bSJed Brown x[0] = tmp[0]; 87*c4762a1bSJed Brown x[1] = tmp[1]; 88*c4762a1bSJed Brown ierr = VecRestoreArray(mu[0],&y);CHKERRQ(ierr); 89*c4762a1bSJed Brown ierr = VecRestoreArray(lambda[0],&x);CHKERRQ(ierr); 90*c4762a1bSJed Brown 91*c4762a1bSJed Brown ierr = VecGetArray(lambda[1],&x);CHKERRQ(ierr); 92*c4762a1bSJed Brown ierr = VecGetArray(mu[1],&y);CHKERRQ(ierr); 93*c4762a1bSJed Brown tmp[0] = A1[0][0]*x[0]+A1[0][1]*x[1]; 94*c4762a1bSJed Brown tmp[1] = A1[1][0]*x[0]+A1[1][1]*x[1]; 95*c4762a1bSJed Brown y[0] = y[0] + A2[0]*x[0]+A2[1]*x[1]; 96*c4762a1bSJed Brown x[0] = tmp[0]; 97*c4762a1bSJed Brown x[1] = tmp[1]; 98*c4762a1bSJed Brown ierr = VecRestoreArray(mu[1],&y);CHKERRQ(ierr); 99*c4762a1bSJed Brown ierr = VecRestoreArray(lambda[1],&x);CHKERRQ(ierr); 100*c4762a1bSJed Brown PetscFunctionReturn(0); 101*c4762a1bSJed Brown } 102*c4762a1bSJed Brown 103*c4762a1bSJed Brown PetscErrorCode PostEventFunction(TS ts,PetscInt nevents,PetscInt event_list[],PetscReal t,Vec U,PetscBool forwardsolve,void* ctx) 104*c4762a1bSJed Brown { 105*c4762a1bSJed Brown AppCtx *actx=(AppCtx*)ctx; 106*c4762a1bSJed Brown PetscErrorCode ierr; 107*c4762a1bSJed Brown 108*c4762a1bSJed Brown PetscFunctionBegin; 109*c4762a1bSJed Brown if (!forwardsolve) { 110*c4762a1bSJed Brown ierr = ShiftGradients(ts,U,actx);CHKERRQ(ierr); 111*c4762a1bSJed Brown } 112*c4762a1bSJed Brown if (actx->mode == 1) { 113*c4762a1bSJed Brown actx->mode = 2; 114*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Change from mode 1 to 2 at t = %f \n",(double)t);CHKERRQ(ierr); 115*c4762a1bSJed Brown } else if (actx->mode == 2) { 116*c4762a1bSJed Brown actx->mode = 1; 117*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Change from mode 2 to 1 at t = %f \n",(double)t);CHKERRQ(ierr); 118*c4762a1bSJed Brown } 119*c4762a1bSJed Brown PetscFunctionReturn(0); 120*c4762a1bSJed Brown } 121*c4762a1bSJed Brown 122*c4762a1bSJed Brown /* 123*c4762a1bSJed Brown Defines the ODE passed to the ODE solver 124*c4762a1bSJed Brown */ 125*c4762a1bSJed Brown static PetscErrorCode IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,void *ctx) 126*c4762a1bSJed Brown { 127*c4762a1bSJed Brown AppCtx *actx=(AppCtx*)ctx; 128*c4762a1bSJed Brown PetscErrorCode ierr; 129*c4762a1bSJed Brown PetscScalar *f; 130*c4762a1bSJed Brown const PetscScalar *u,*udot; 131*c4762a1bSJed Brown 132*c4762a1bSJed Brown PetscFunctionBegin; 133*c4762a1bSJed Brown /* The next three lines allow us to access the entries of the vectors directly */ 134*c4762a1bSJed Brown ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 135*c4762a1bSJed Brown ierr = VecGetArrayRead(Udot,&udot);CHKERRQ(ierr); 136*c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 137*c4762a1bSJed Brown 138*c4762a1bSJed Brown if (actx->mode == 1) { 139*c4762a1bSJed Brown f[0] = udot[0]-u[0]+100*u[1]; 140*c4762a1bSJed Brown f[1] = udot[1]-10*u[0]-u[1]; 141*c4762a1bSJed Brown } else if (actx->mode == 2) { 142*c4762a1bSJed Brown f[0] = udot[0]-u[0]-10*u[1]; 143*c4762a1bSJed Brown f[1] = udot[1]+100*u[0]-u[1]; 144*c4762a1bSJed Brown } 145*c4762a1bSJed Brown 146*c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 147*c4762a1bSJed Brown ierr = VecRestoreArrayRead(Udot,&udot);CHKERRQ(ierr); 148*c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 149*c4762a1bSJed Brown PetscFunctionReturn(0); 150*c4762a1bSJed Brown } 151*c4762a1bSJed Brown 152*c4762a1bSJed Brown /* 153*c4762a1bSJed Brown Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian. 154*c4762a1bSJed Brown */ 155*c4762a1bSJed Brown static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,void *ctx) 156*c4762a1bSJed Brown { 157*c4762a1bSJed Brown AppCtx *actx=(AppCtx*)ctx; 158*c4762a1bSJed Brown PetscErrorCode ierr; 159*c4762a1bSJed Brown PetscInt rowcol[] = {0,1}; 160*c4762a1bSJed Brown PetscScalar J[2][2]; 161*c4762a1bSJed Brown const PetscScalar *u,*udot; 162*c4762a1bSJed Brown 163*c4762a1bSJed Brown PetscFunctionBegin; 164*c4762a1bSJed Brown ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 165*c4762a1bSJed Brown ierr = VecGetArrayRead(Udot,&udot);CHKERRQ(ierr); 166*c4762a1bSJed Brown 167*c4762a1bSJed Brown if (actx->mode == 1) { 168*c4762a1bSJed Brown J[0][0] = a-1; J[0][1] = 100; 169*c4762a1bSJed Brown J[1][0] = -10; J[1][1] = a-1; 170*c4762a1bSJed Brown } else if (actx->mode == 2) { 171*c4762a1bSJed Brown J[0][0] = a-1; J[0][1] = -10; 172*c4762a1bSJed Brown J[1][0] = 100; J[1][1] = a-1; 173*c4762a1bSJed Brown } 174*c4762a1bSJed Brown ierr = MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);CHKERRQ(ierr); 175*c4762a1bSJed Brown 176*c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 177*c4762a1bSJed Brown ierr = VecRestoreArrayRead(Udot,&udot);CHKERRQ(ierr); 178*c4762a1bSJed Brown 179*c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 180*c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 181*c4762a1bSJed Brown if (A != B) { 182*c4762a1bSJed Brown ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 183*c4762a1bSJed Brown ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 184*c4762a1bSJed Brown } 185*c4762a1bSJed Brown PetscFunctionReturn(0); 186*c4762a1bSJed Brown } 187*c4762a1bSJed Brown 188*c4762a1bSJed Brown int main(int argc,char **argv) 189*c4762a1bSJed Brown { 190*c4762a1bSJed Brown TS ts; /* ODE integrator */ 191*c4762a1bSJed Brown Vec U; /* solution will be stored here */ 192*c4762a1bSJed Brown Mat A; /* Jacobian matrix */ 193*c4762a1bSJed Brown Mat Ap; /* dfdp */ 194*c4762a1bSJed Brown PetscErrorCode ierr; 195*c4762a1bSJed Brown PetscMPIInt size; 196*c4762a1bSJed Brown PetscInt n = 2; 197*c4762a1bSJed Brown PetscScalar *u; 198*c4762a1bSJed Brown AppCtx app; 199*c4762a1bSJed Brown PetscInt direction[1]; 200*c4762a1bSJed Brown PetscBool terminate[1]; 201*c4762a1bSJed Brown PetscReal delta,tmp[2],sensi[2]; 202*c4762a1bSJed Brown 203*c4762a1bSJed Brown delta = 1e-8; 204*c4762a1bSJed Brown 205*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 206*c4762a1bSJed Brown Initialize program 207*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 208*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 209*c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 210*c4762a1bSJed Brown if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs"); 211*c4762a1bSJed Brown app.mode = 1; 212*c4762a1bSJed Brown app.lambda1 = 2.75; 213*c4762a1bSJed Brown app.lambda2 = 0.36; 214*c4762a1bSJed Brown ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"ex1 options","");CHKERRQ(ierr); 215*c4762a1bSJed Brown { 216*c4762a1bSJed Brown ierr = PetscOptionsReal("-lambda1","","",app.lambda1,&app.lambda1,NULL);CHKERRQ(ierr); 217*c4762a1bSJed Brown ierr = PetscOptionsReal("-lambda2","","",app.lambda2,&app.lambda2,NULL);CHKERRQ(ierr); 218*c4762a1bSJed Brown } 219*c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 220*c4762a1bSJed Brown 221*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 222*c4762a1bSJed Brown Create necessary matrix and vectors 223*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 224*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 225*c4762a1bSJed Brown ierr = MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 226*c4762a1bSJed Brown ierr = MatSetType(A,MATDENSE);CHKERRQ(ierr); 227*c4762a1bSJed Brown ierr = MatSetFromOptions(A);CHKERRQ(ierr); 228*c4762a1bSJed Brown ierr = MatSetUp(A);CHKERRQ(ierr); 229*c4762a1bSJed Brown 230*c4762a1bSJed Brown ierr = MatCreateVecs(A,&U,NULL);CHKERRQ(ierr); 231*c4762a1bSJed Brown 232*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&Ap);CHKERRQ(ierr); 233*c4762a1bSJed Brown ierr = MatSetSizes(Ap,n,1,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 234*c4762a1bSJed Brown ierr = MatSetType(Ap,MATDENSE);CHKERRQ(ierr); 235*c4762a1bSJed Brown ierr = MatSetFromOptions(Ap);CHKERRQ(ierr); 236*c4762a1bSJed Brown ierr = MatSetUp(Ap);CHKERRQ(ierr); 237*c4762a1bSJed Brown ierr = MatZeroEntries(Ap);CHKERRQ(ierr); /* initialize to zeros */ 238*c4762a1bSJed Brown 239*c4762a1bSJed Brown ierr = VecGetArray(U,&u);CHKERRQ(ierr); 240*c4762a1bSJed Brown u[0] = 0; 241*c4762a1bSJed Brown u[1] = 1; 242*c4762a1bSJed Brown ierr = VecRestoreArray(U,&u);CHKERRQ(ierr); 243*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 244*c4762a1bSJed Brown Create timestepping solver context 245*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 246*c4762a1bSJed Brown ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr); 247*c4762a1bSJed Brown ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr); 248*c4762a1bSJed Brown ierr = TSSetType(ts,TSCN);CHKERRQ(ierr); 249*c4762a1bSJed Brown ierr = TSSetIFunction(ts,NULL,(TSIFunction)IFunction,&app);CHKERRQ(ierr); 250*c4762a1bSJed Brown ierr = TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,&app);CHKERRQ(ierr); 251*c4762a1bSJed Brown 252*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 253*c4762a1bSJed Brown Set initial conditions 254*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 255*c4762a1bSJed Brown ierr = TSSetSolution(ts,U);CHKERRQ(ierr); 256*c4762a1bSJed Brown 257*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 258*c4762a1bSJed Brown Set solver options 259*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 260*c4762a1bSJed Brown ierr = TSSetMaxTime(ts,0.125);CHKERRQ(ierr); 261*c4762a1bSJed Brown ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr); 262*c4762a1bSJed Brown ierr = TSSetTimeStep(ts,1./256.);CHKERRQ(ierr); 263*c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 264*c4762a1bSJed Brown 265*c4762a1bSJed Brown /* Set directions and terminate flags for the two events */ 266*c4762a1bSJed Brown direction[0] = 0; 267*c4762a1bSJed Brown terminate[0] = PETSC_FALSE; 268*c4762a1bSJed Brown ierr = TSSetEventHandler(ts,1,direction,terminate,EventFunction,PostEventFunction,(void*)&app);CHKERRQ(ierr); 269*c4762a1bSJed Brown 270*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 271*c4762a1bSJed Brown Run timestepping solver 272*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 273*c4762a1bSJed Brown ierr = TSSolve(ts,U);CHKERRQ(ierr); 274*c4762a1bSJed Brown 275*c4762a1bSJed Brown ierr = VecGetArray(U,&u);CHKERRQ(ierr); 276*c4762a1bSJed Brown tmp[0] = u[0]; 277*c4762a1bSJed Brown tmp[1] = u[1]; 278*c4762a1bSJed Brown 279*c4762a1bSJed Brown u[0] = 0+delta; 280*c4762a1bSJed Brown u[1] = 1; 281*c4762a1bSJed Brown ierr = VecRestoreArray(U,&u);CHKERRQ(ierr); 282*c4762a1bSJed Brown 283*c4762a1bSJed Brown ierr = FWDRun(ts,U,(void*)&app);CHKERRQ(ierr); 284*c4762a1bSJed Brown 285*c4762a1bSJed Brown ierr = VecGetArray(U,&u);CHKERRQ(ierr); 286*c4762a1bSJed Brown sensi[0] = (u[0]-tmp[0])/delta; 287*c4762a1bSJed Brown sensi[1] = (u[1]-tmp[1])/delta; 288*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"d x1(tf) /d x1(t0) = %f d x2(tf) / d x1(t0) = %f \n",sensi[0],sensi[1]);CHKERRQ(ierr); 289*c4762a1bSJed Brown u[0] = 0; 290*c4762a1bSJed Brown u[1] = 1+delta; 291*c4762a1bSJed Brown ierr = VecRestoreArray(U,&u);CHKERRQ(ierr); 292*c4762a1bSJed Brown 293*c4762a1bSJed Brown ierr = FWDRun(ts,U,(void*)&app);CHKERRQ(ierr); 294*c4762a1bSJed Brown 295*c4762a1bSJed Brown ierr = VecGetArray(U,&u);CHKERRQ(ierr); 296*c4762a1bSJed Brown sensi[0] = (u[0]-tmp[0])/delta; 297*c4762a1bSJed Brown sensi[1] = (u[1]-tmp[1])/delta; 298*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"d x1(tf) /d x2(t0) = %f d x2(tf) / d x2(t0) = %f \n",sensi[0],sensi[1]);CHKERRQ(ierr); 299*c4762a1bSJed Brown u[0] = 0; 300*c4762a1bSJed Brown u[1] = 1; 301*c4762a1bSJed Brown app.lambda1 = app.lambda1+delta; 302*c4762a1bSJed Brown ierr = VecRestoreArray(U,&u);CHKERRQ(ierr); 303*c4762a1bSJed Brown 304*c4762a1bSJed Brown ierr = FWDRun(ts,U,(void*)&app);CHKERRQ(ierr); 305*c4762a1bSJed Brown 306*c4762a1bSJed Brown ierr = VecGetArray(U,&u);CHKERRQ(ierr); 307*c4762a1bSJed Brown sensi[0] = (u[0]-tmp[0])/delta; 308*c4762a1bSJed Brown sensi[1] = (u[1]-tmp[1])/delta; 309*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Final gradients: d x1(tf) /d p = %f d x2(tf) / d p = %f \n",sensi[0],sensi[1]);CHKERRQ(ierr); 310*c4762a1bSJed Brown ierr = VecRestoreArray(U,&u);CHKERRQ(ierr); 311*c4762a1bSJed Brown 312*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 313*c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they are no longer needed. 314*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 315*c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 316*c4762a1bSJed Brown ierr = VecDestroy(&U);CHKERRQ(ierr); 317*c4762a1bSJed Brown ierr = TSDestroy(&ts);CHKERRQ(ierr); 318*c4762a1bSJed Brown 319*c4762a1bSJed Brown ierr = MatDestroy(&Ap);CHKERRQ(ierr); 320*c4762a1bSJed Brown ierr = PetscFinalize(); 321*c4762a1bSJed Brown return ierr; 322*c4762a1bSJed Brown } 323*c4762a1bSJed Brown 324*c4762a1bSJed Brown PetscErrorCode FWDRun(TS ts, Vec U0, void *ctx0) 325*c4762a1bSJed Brown { 326*c4762a1bSJed Brown Vec U; /* solution will be stored here */ 327*c4762a1bSJed Brown PetscErrorCode ierr; 328*c4762a1bSJed Brown AppCtx *ctx=(AppCtx*)ctx0; 329*c4762a1bSJed Brown 330*c4762a1bSJed Brown PetscFunctionBeginUser; 331*c4762a1bSJed Brown ierr = TSGetSolution(ts,&U);CHKERRQ(ierr); 332*c4762a1bSJed Brown ierr = VecCopy(U0,U);CHKERRQ(ierr); 333*c4762a1bSJed Brown 334*c4762a1bSJed Brown ctx->mode = 1; 335*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 336*c4762a1bSJed Brown Run timestepping solver 337*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 338*c4762a1bSJed Brown ierr = TSSetTime(ts, 0.0);CHKERRQ(ierr); 339*c4762a1bSJed Brown 340*c4762a1bSJed Brown ierr = TSSolve(ts,U);CHKERRQ(ierr); 341*c4762a1bSJed Brown 342*c4762a1bSJed Brown PetscFunctionReturn(0); 343*c4762a1bSJed Brown } 344*c4762a1bSJed Brown 345*c4762a1bSJed Brown 346*c4762a1bSJed Brown /*TEST 347*c4762a1bSJed Brown 348*c4762a1bSJed Brown build: 349*c4762a1bSJed Brown requires: !define(PETSC_USE_CXXCOMPLEX) 350*c4762a1bSJed Brown 351*c4762a1bSJed Brown test: 352*c4762a1bSJed Brown args: -ts_event_tol 1e-9 353*c4762a1bSJed Brown timeoutfactor: 18 354*c4762a1bSJed Brown requires: !single 355*c4762a1bSJed Brown 356*c4762a1bSJed Brown TEST*/ 357