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