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