1*c4762a1bSJed Brown static char help[] = "Performs adjoint sensitivity analysis for the van der Pol equation.\n"; 2*c4762a1bSJed Brown 3*c4762a1bSJed Brown /* 4*c4762a1bSJed Brown Concepts: TS^time-dependent nonlinear problems 5*c4762a1bSJed Brown Concepts: TS^van der Pol equation DAE equivalent 6*c4762a1bSJed Brown Concepts: TS^adjoint sensitivity analysis 7*c4762a1bSJed Brown Processors: 1 8*c4762a1bSJed Brown */ 9*c4762a1bSJed Brown /* ------------------------------------------------------------------------ 10*c4762a1bSJed Brown 11*c4762a1bSJed Brown This program solves the van der Pol DAE ODE equivalent 12*c4762a1bSJed Brown [ u_1' ] = [ u_2 ] (2) 13*c4762a1bSJed Brown [ u_2' ] [ \mu ((1 - u_1^2) u_2 - u_1) ] 14*c4762a1bSJed Brown on the domain 0 <= x <= 1, with the boundary conditions 15*c4762a1bSJed Brown u_1(0) = 2, u_2(0) = - 2/3 +10/(81*\mu) - 292/(2187*\mu^2), 16*c4762a1bSJed Brown and 17*c4762a1bSJed Brown \mu = 10^6 ( y'(0) ~ -0.6666665432100101)., 18*c4762a1bSJed Brown and computes the sensitivities of the final solution w.r.t. initial conditions and parameter \mu with the implicit theta method and its discrete adjoint. 19*c4762a1bSJed Brown 20*c4762a1bSJed Brown Notes: 21*c4762a1bSJed Brown This code demonstrates the TSAdjoint interface to a DAE system. 22*c4762a1bSJed Brown 23*c4762a1bSJed Brown The user provides the implicit right-hand-side function 24*c4762a1bSJed Brown [ G(u',u,t) ] = [u' - f(u,t)] = [ u_1'] - [ u_2 ] 25*c4762a1bSJed Brown [ u_2'] [ \mu ((1-u_1^2)u_2-u_1) ] 26*c4762a1bSJed Brown 27*c4762a1bSJed Brown and the Jacobian of G (from the PETSc user manual) 28*c4762a1bSJed Brown 29*c4762a1bSJed Brown dG dG 30*c4762a1bSJed Brown J(G) = a * -- + -- 31*c4762a1bSJed Brown du' du 32*c4762a1bSJed Brown 33*c4762a1bSJed Brown and the JacobianP of the explicit right-hand side of (2) f(u,t) ( which is equivalent to -G(0,u,t) ). 34*c4762a1bSJed Brown df [ 0 ] 35*c4762a1bSJed Brown -- = [ ] 36*c4762a1bSJed Brown dp [ (1 - u_1^2) u_2 - u_1 ]. 37*c4762a1bSJed Brown 38*c4762a1bSJed Brown See ex20.c for more details on the Jacobian. 39*c4762a1bSJed Brown 40*c4762a1bSJed Brown ------------------------------------------------------------------------- */ 41*c4762a1bSJed Brown #include <petscts.h> 42*c4762a1bSJed Brown #include <petsctao.h> 43*c4762a1bSJed Brown 44*c4762a1bSJed Brown typedef struct _n_User *User; 45*c4762a1bSJed Brown struct _n_User { 46*c4762a1bSJed Brown PetscReal mu; 47*c4762a1bSJed Brown PetscReal next_output; 48*c4762a1bSJed Brown 49*c4762a1bSJed Brown /* Sensitivity analysis support */ 50*c4762a1bSJed Brown PetscInt steps; 51*c4762a1bSJed Brown PetscReal ftime; 52*c4762a1bSJed Brown Mat A; /* Jacobian matrix */ 53*c4762a1bSJed Brown Mat Jacp; /* JacobianP matrix */ 54*c4762a1bSJed Brown Vec U,lambda[2],mup[2]; /* adjoint variables */ 55*c4762a1bSJed Brown }; 56*c4762a1bSJed Brown 57*c4762a1bSJed Brown /* ----------------------- Explicit form of the ODE -------------------- */ 58*c4762a1bSJed Brown 59*c4762a1bSJed Brown static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec U,Vec F,void *ctx) 60*c4762a1bSJed Brown { 61*c4762a1bSJed Brown PetscErrorCode ierr; 62*c4762a1bSJed Brown User user = (User)ctx; 63*c4762a1bSJed Brown PetscScalar *f; 64*c4762a1bSJed Brown const PetscScalar *u; 65*c4762a1bSJed Brown 66*c4762a1bSJed Brown PetscFunctionBeginUser; 67*c4762a1bSJed Brown ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 68*c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 69*c4762a1bSJed Brown f[0] = u[1]; 70*c4762a1bSJed Brown f[1] = user->mu*((1.-u[0]*u[0])*u[1]-u[0]); 71*c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 72*c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 73*c4762a1bSJed Brown PetscFunctionReturn(0); 74*c4762a1bSJed Brown } 75*c4762a1bSJed Brown 76*c4762a1bSJed Brown static PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec U,Mat A,Mat B,void *ctx) 77*c4762a1bSJed Brown { 78*c4762a1bSJed Brown PetscErrorCode ierr; 79*c4762a1bSJed Brown User user = (User)ctx; 80*c4762a1bSJed Brown PetscReal mu = user->mu; 81*c4762a1bSJed Brown PetscInt rowcol[] = {0,1}; 82*c4762a1bSJed Brown PetscScalar J[2][2]; 83*c4762a1bSJed Brown const PetscScalar *u; 84*c4762a1bSJed Brown 85*c4762a1bSJed Brown PetscFunctionBeginUser; 86*c4762a1bSJed Brown ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 87*c4762a1bSJed Brown J[0][0] = 0; 88*c4762a1bSJed Brown J[1][0] = -mu*(2.0*u[1]*u[0]+1.); 89*c4762a1bSJed Brown J[0][1] = 1.0; 90*c4762a1bSJed Brown J[1][1] = mu*(1.0-u[0]*u[0]); 91*c4762a1bSJed Brown ierr = MatSetValues(A,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);CHKERRQ(ierr); 92*c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 93*c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 94*c4762a1bSJed Brown if (A != B) { 95*c4762a1bSJed Brown ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 96*c4762a1bSJed Brown ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 97*c4762a1bSJed Brown } 98*c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 99*c4762a1bSJed Brown PetscFunctionReturn(0); 100*c4762a1bSJed Brown } 101*c4762a1bSJed Brown 102*c4762a1bSJed Brown /* ----------------------- Implicit form of the ODE -------------------- */ 103*c4762a1bSJed Brown 104*c4762a1bSJed Brown static PetscErrorCode IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,void *ctx) 105*c4762a1bSJed Brown { 106*c4762a1bSJed Brown PetscErrorCode ierr; 107*c4762a1bSJed Brown User user = (User)ctx; 108*c4762a1bSJed Brown const PetscScalar *u,*udot; 109*c4762a1bSJed Brown PetscScalar *f; 110*c4762a1bSJed Brown 111*c4762a1bSJed Brown PetscFunctionBeginUser; 112*c4762a1bSJed Brown ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 113*c4762a1bSJed Brown ierr = VecGetArrayRead(Udot,&udot);CHKERRQ(ierr); 114*c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 115*c4762a1bSJed Brown f[0] = udot[0] - u[1]; 116*c4762a1bSJed Brown f[1] = udot[1] - user->mu*((1.0-u[0]*u[0])*u[1] - u[0]); 117*c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 118*c4762a1bSJed Brown ierr = VecRestoreArrayRead(Udot,&udot);CHKERRQ(ierr); 119*c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 120*c4762a1bSJed Brown PetscFunctionReturn(0); 121*c4762a1bSJed Brown } 122*c4762a1bSJed Brown 123*c4762a1bSJed Brown static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,void *ctx) 124*c4762a1bSJed Brown { 125*c4762a1bSJed Brown PetscErrorCode ierr; 126*c4762a1bSJed Brown User user = (User)ctx; 127*c4762a1bSJed Brown PetscInt rowcol[] = {0,1}; 128*c4762a1bSJed Brown PetscScalar J[2][2]; 129*c4762a1bSJed Brown const PetscScalar *u; 130*c4762a1bSJed Brown 131*c4762a1bSJed Brown PetscFunctionBeginUser; 132*c4762a1bSJed Brown ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 133*c4762a1bSJed Brown 134*c4762a1bSJed Brown J[0][0] = a; J[0][1] = -1.0; 135*c4762a1bSJed Brown J[1][0] = user->mu*(2.0*u[0]*u[1] + 1.0); J[1][1] = a - user->mu*(1.0-u[0]*u[0]); 136*c4762a1bSJed Brown 137*c4762a1bSJed Brown ierr = MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);CHKERRQ(ierr); 138*c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 139*c4762a1bSJed Brown 140*c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 141*c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 142*c4762a1bSJed Brown if (B && A != B) { 143*c4762a1bSJed Brown ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 144*c4762a1bSJed Brown ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 145*c4762a1bSJed Brown } 146*c4762a1bSJed Brown PetscFunctionReturn(0); 147*c4762a1bSJed Brown } 148*c4762a1bSJed Brown 149*c4762a1bSJed Brown static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec U,Mat A,void *ctx) 150*c4762a1bSJed Brown { 151*c4762a1bSJed Brown PetscErrorCode ierr; 152*c4762a1bSJed Brown PetscInt row[] = {0,1},col[]={0}; 153*c4762a1bSJed Brown PetscScalar J[2][1]; 154*c4762a1bSJed Brown const PetscScalar *u; 155*c4762a1bSJed Brown 156*c4762a1bSJed Brown PetscFunctionBeginUser; 157*c4762a1bSJed Brown ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 158*c4762a1bSJed Brown J[0][0] = 0; 159*c4762a1bSJed Brown J[1][0] = (1.-u[0]*u[0])*u[1]-u[0]; 160*c4762a1bSJed Brown ierr = MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES);CHKERRQ(ierr); 161*c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 162*c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 163*c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 164*c4762a1bSJed Brown PetscFunctionReturn(0); 165*c4762a1bSJed Brown } 166*c4762a1bSJed Brown 167*c4762a1bSJed Brown /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */ 168*c4762a1bSJed Brown static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec U,void *ctx) 169*c4762a1bSJed Brown { 170*c4762a1bSJed Brown PetscErrorCode ierr; 171*c4762a1bSJed Brown const PetscScalar *u; 172*c4762a1bSJed Brown PetscReal tfinal, dt; 173*c4762a1bSJed Brown User user = (User)ctx; 174*c4762a1bSJed Brown Vec interpolatedU; 175*c4762a1bSJed Brown 176*c4762a1bSJed Brown PetscFunctionBeginUser; 177*c4762a1bSJed Brown ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr); 178*c4762a1bSJed Brown ierr = TSGetMaxTime(ts,&tfinal);CHKERRQ(ierr); 179*c4762a1bSJed Brown 180*c4762a1bSJed Brown while (user->next_output <= t && user->next_output <= tfinal) { 181*c4762a1bSJed Brown ierr = VecDuplicate(U,&interpolatedU);CHKERRQ(ierr); 182*c4762a1bSJed Brown ierr = TSInterpolate(ts,user->next_output,interpolatedU);CHKERRQ(ierr); 183*c4762a1bSJed Brown ierr = VecGetArrayRead(interpolatedU,&u);CHKERRQ(ierr); 184*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"[%g] %D TS %g (dt = %g) X %g %g\n", 185*c4762a1bSJed Brown (double)user->next_output,step,(double)t,(double)dt,(double)PetscRealPart(u[0]), 186*c4762a1bSJed Brown (double)PetscRealPart(u[1]));CHKERRQ(ierr); 187*c4762a1bSJed Brown ierr = VecRestoreArrayRead(interpolatedU,&u);CHKERRQ(ierr); 188*c4762a1bSJed Brown ierr = VecDestroy(&interpolatedU);CHKERRQ(ierr); 189*c4762a1bSJed Brown user->next_output += 0.1; 190*c4762a1bSJed Brown } 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; 197*c4762a1bSJed Brown PetscBool monitor = PETSC_FALSE,implicitform = PETSC_TRUE; 198*c4762a1bSJed Brown PetscScalar *x_ptr,*y_ptr,derp; 199*c4762a1bSJed Brown PetscMPIInt size; 200*c4762a1bSJed Brown struct _n_User user; 201*c4762a1bSJed Brown PetscErrorCode ierr; 202*c4762a1bSJed Brown 203*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 204*c4762a1bSJed Brown Initialize program 205*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 206*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr; 207*c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 208*c4762a1bSJed Brown if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 209*c4762a1bSJed Brown 210*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 211*c4762a1bSJed Brown Set runtime options 212*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 213*c4762a1bSJed Brown user.next_output = 0.0; 214*c4762a1bSJed Brown user.mu = 1.0e3; 215*c4762a1bSJed Brown user.steps = 0; 216*c4762a1bSJed Brown user.ftime = 0.5; 217*c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL);CHKERRQ(ierr); 218*c4762a1bSJed Brown ierr = PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,NULL);CHKERRQ(ierr); 219*c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-implicitform",&implicitform,NULL);CHKERRQ(ierr); 220*c4762a1bSJed Brown 221*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 222*c4762a1bSJed Brown Create necessary matrix and vectors, solve same ODE on every process 223*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 224*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&user.A);CHKERRQ(ierr); 225*c4762a1bSJed Brown ierr = MatSetSizes(user.A,PETSC_DECIDE,PETSC_DECIDE,2,2);CHKERRQ(ierr); 226*c4762a1bSJed Brown ierr = MatSetFromOptions(user.A);CHKERRQ(ierr); 227*c4762a1bSJed Brown ierr = MatSetUp(user.A);CHKERRQ(ierr); 228*c4762a1bSJed Brown ierr = MatCreateVecs(user.A,&user.U,NULL);CHKERRQ(ierr); 229*c4762a1bSJed Brown 230*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&user.Jacp);CHKERRQ(ierr); 231*c4762a1bSJed Brown ierr = MatSetSizes(user.Jacp,PETSC_DECIDE,PETSC_DECIDE,2,1);CHKERRQ(ierr); 232*c4762a1bSJed Brown ierr = MatSetFromOptions(user.Jacp);CHKERRQ(ierr); 233*c4762a1bSJed Brown ierr = MatSetUp(user.Jacp);CHKERRQ(ierr); 234*c4762a1bSJed Brown 235*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 236*c4762a1bSJed Brown Create timestepping solver context 237*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 238*c4762a1bSJed Brown ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr); 239*c4762a1bSJed Brown ierr = TSSetEquationType(ts,TS_EQ_ODE_EXPLICIT);CHKERRQ(ierr); /* less Jacobian evaluations when adjoint BEuler is used, otherwise no effect */ 240*c4762a1bSJed Brown if (implicitform) { 241*c4762a1bSJed Brown ierr = TSSetIFunction(ts,NULL,IFunction,&user);CHKERRQ(ierr); 242*c4762a1bSJed Brown ierr = TSSetIJacobian(ts,user.A,user.A,IJacobian,&user);CHKERRQ(ierr); 243*c4762a1bSJed Brown ierr = TSSetType(ts,TSCN);CHKERRQ(ierr); 244*c4762a1bSJed Brown } else { 245*c4762a1bSJed Brown ierr = TSSetRHSFunction(ts,NULL,RHSFunction,&user);CHKERRQ(ierr); 246*c4762a1bSJed Brown ierr = TSSetRHSJacobian(ts,user.A,user.A,RHSJacobian,&user);CHKERRQ(ierr); 247*c4762a1bSJed Brown ierr = TSSetType(ts,TSRK);CHKERRQ(ierr); 248*c4762a1bSJed Brown } 249*c4762a1bSJed Brown ierr = TSSetRHSJacobianP(ts,user.Jacp,RHSJacobianP,&user);CHKERRQ(ierr); 250*c4762a1bSJed Brown ierr = TSSetMaxTime(ts,user.ftime);CHKERRQ(ierr); 251*c4762a1bSJed Brown ierr = TSSetTimeStep(ts,0.001);CHKERRQ(ierr); 252*c4762a1bSJed Brown ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr); 253*c4762a1bSJed Brown if (monitor) { 254*c4762a1bSJed Brown ierr = TSMonitorSet(ts,Monitor,&user,NULL);CHKERRQ(ierr); 255*c4762a1bSJed Brown } 256*c4762a1bSJed Brown 257*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 258*c4762a1bSJed Brown Set initial conditions 259*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 260*c4762a1bSJed Brown ierr = VecGetArray(user.U,&x_ptr);CHKERRQ(ierr); 261*c4762a1bSJed Brown x_ptr[0] = 2.0; 262*c4762a1bSJed Brown x_ptr[1] = -2.0/3.0 + 10.0/(81.0*user.mu) - 292.0/(2187.0*user.mu*user.mu); 263*c4762a1bSJed Brown ierr = VecRestoreArray(user.U,&x_ptr);CHKERRQ(ierr); 264*c4762a1bSJed Brown ierr = TSSetTimeStep(ts,0.001);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 runtime options 273*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 274*c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 275*c4762a1bSJed Brown 276*c4762a1bSJed Brown ierr = TSSolve(ts,user.U);CHKERRQ(ierr); 277*c4762a1bSJed Brown ierr = TSGetSolveTime(ts,&user.ftime);CHKERRQ(ierr); 278*c4762a1bSJed Brown ierr = TSGetStepNumber(ts,&user.steps);CHKERRQ(ierr); 279*c4762a1bSJed Brown 280*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 281*c4762a1bSJed Brown Adjoint model starts here 282*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 283*c4762a1bSJed Brown ierr = MatCreateVecs(user.A,&user.lambda[0],NULL);CHKERRQ(ierr); 284*c4762a1bSJed Brown /* Set initial conditions for the adjoint integration */ 285*c4762a1bSJed Brown ierr = VecGetArray(user.lambda[0],&y_ptr);CHKERRQ(ierr); 286*c4762a1bSJed Brown y_ptr[0] = 1.0; y_ptr[1] = 0.0; 287*c4762a1bSJed Brown ierr = VecRestoreArray(user.lambda[0],&y_ptr);CHKERRQ(ierr); 288*c4762a1bSJed Brown ierr = MatCreateVecs(user.A,&user.lambda[1],NULL);CHKERRQ(ierr); 289*c4762a1bSJed Brown ierr = VecGetArray(user.lambda[1],&y_ptr);CHKERRQ(ierr); 290*c4762a1bSJed Brown y_ptr[0] = 0.0; y_ptr[1] = 1.0; 291*c4762a1bSJed Brown ierr = VecRestoreArray(user.lambda[1],&y_ptr);CHKERRQ(ierr); 292*c4762a1bSJed Brown 293*c4762a1bSJed Brown ierr = MatCreateVecs(user.Jacp,&user.mup[0],NULL);CHKERRQ(ierr); 294*c4762a1bSJed Brown ierr = VecGetArray(user.mup[0],&x_ptr);CHKERRQ(ierr); 295*c4762a1bSJed Brown x_ptr[0] = 0.0; 296*c4762a1bSJed Brown ierr = VecRestoreArray(user.mup[0],&x_ptr);CHKERRQ(ierr); 297*c4762a1bSJed Brown ierr = MatCreateVecs(user.Jacp,&user.mup[1],NULL);CHKERRQ(ierr); 298*c4762a1bSJed Brown ierr = VecGetArray(user.mup[1],&x_ptr);CHKERRQ(ierr); 299*c4762a1bSJed Brown x_ptr[0] = 0.0; 300*c4762a1bSJed Brown ierr = VecRestoreArray(user.mup[1],&x_ptr);CHKERRQ(ierr); 301*c4762a1bSJed Brown 302*c4762a1bSJed Brown ierr = TSSetCostGradients(ts,2,user.lambda,user.mup);CHKERRQ(ierr); 303*c4762a1bSJed Brown 304*c4762a1bSJed Brown ierr = TSAdjointSolve(ts);CHKERRQ(ierr); 305*c4762a1bSJed Brown 306*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt initial conditions: d[y(tf)]/d[y0] d[y(tf)]/d[z0]\n");CHKERRQ(ierr); 307*c4762a1bSJed Brown ierr = VecView(user.lambda[0],PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 308*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt initial conditions: d[z(tf)]/d[y0] d[z(tf)]/d[z0]\n");CHKERRQ(ierr); 309*c4762a1bSJed Brown ierr = VecView(user.lambda[1],PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 310*c4762a1bSJed Brown 311*c4762a1bSJed Brown ierr = VecGetArray(user.mup[0],&x_ptr);CHKERRQ(ierr); 312*c4762a1bSJed Brown ierr = VecGetArray(user.lambda[0],&y_ptr);CHKERRQ(ierr); 313*c4762a1bSJed Brown derp = y_ptr[1]*(-10.0/(81.0*user.mu*user.mu)+2.0*292.0/(2187.0*user.mu*user.mu*user.mu))+x_ptr[0]; 314*c4762a1bSJed Brown ierr = VecRestoreArray(user.mup[0],&x_ptr);CHKERRQ(ierr); 315*c4762a1bSJed Brown ierr = VecRestoreArray(user.lambda[0],&y_ptr);CHKERRQ(ierr); 316*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt parameters: d[y(tf)]/d[mu]\n%g\n",(double)PetscRealPart(derp));CHKERRQ(ierr); 317*c4762a1bSJed Brown 318*c4762a1bSJed Brown ierr = VecGetArray(user.mup[1],&x_ptr);CHKERRQ(ierr); 319*c4762a1bSJed Brown ierr = VecGetArray(user.lambda[1],&y_ptr);CHKERRQ(ierr); 320*c4762a1bSJed Brown derp = y_ptr[1]*(-10.0/(81.0*user.mu*user.mu)+2.0*292.0/(2187.0*user.mu*user.mu*user.mu))+x_ptr[0]; 321*c4762a1bSJed Brown ierr = VecRestoreArray(user.mup[1],&x_ptr);CHKERRQ(ierr); 322*c4762a1bSJed Brown ierr = VecRestoreArray(user.lambda[1],&y_ptr);CHKERRQ(ierr); 323*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"\n sensivitity wrt parameters: d[z(tf)]/d[mu]\n%g\n",(double)PetscRealPart(derp));CHKERRQ(ierr); 324*c4762a1bSJed Brown 325*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 326*c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 327*c4762a1bSJed Brown are no longer needed. 328*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 329*c4762a1bSJed Brown ierr = MatDestroy(&user.A);CHKERRQ(ierr); 330*c4762a1bSJed Brown ierr = MatDestroy(&user.Jacp);CHKERRQ(ierr); 331*c4762a1bSJed Brown ierr = VecDestroy(&user.U);CHKERRQ(ierr); 332*c4762a1bSJed Brown ierr = VecDestroy(&user.lambda[0]);CHKERRQ(ierr); 333*c4762a1bSJed Brown ierr = VecDestroy(&user.lambda[1]);CHKERRQ(ierr); 334*c4762a1bSJed Brown ierr = VecDestroy(&user.mup[0]);CHKERRQ(ierr); 335*c4762a1bSJed Brown ierr = VecDestroy(&user.mup[1]);CHKERRQ(ierr); 336*c4762a1bSJed Brown ierr = TSDestroy(&ts);CHKERRQ(ierr); 337*c4762a1bSJed Brown 338*c4762a1bSJed Brown ierr = PetscFinalize(); 339*c4762a1bSJed Brown return(ierr); 340*c4762a1bSJed Brown } 341*c4762a1bSJed Brown 342*c4762a1bSJed Brown /*TEST 343*c4762a1bSJed Brown 344*c4762a1bSJed Brown test: 345*c4762a1bSJed Brown requires: revolve 346*c4762a1bSJed Brown args: -monitor 0 -ts_type theta -ts_theta_endpoint -ts_theta_theta 0.5 -viewer_binary_skip_info -ts_dt 0.001 -mu 100000 347*c4762a1bSJed Brown 348*c4762a1bSJed Brown test: 349*c4762a1bSJed Brown suffix: 2 350*c4762a1bSJed Brown args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_solution_only 351*c4762a1bSJed Brown 352*c4762a1bSJed Brown test: 353*c4762a1bSJed Brown suffix: 3 354*c4762a1bSJed Brown args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_solution_only 0 355*c4762a1bSJed Brown output_file: output/ex20adj_2.out 356*c4762a1bSJed Brown 357*c4762a1bSJed Brown test: 358*c4762a1bSJed Brown suffix: 4 359*c4762a1bSJed Brown args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_stride 5 -ts_trajectory_solution_only -ts_trajectory_save_stack 360*c4762a1bSJed Brown output_file: output/ex20adj_2.out 361*c4762a1bSJed Brown 362*c4762a1bSJed Brown test: 363*c4762a1bSJed Brown suffix: 5 364*c4762a1bSJed Brown args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_stride 5 -ts_trajectory_solution_only 0 -ts_trajectory_save_stack 365*c4762a1bSJed Brown output_file: output/ex20adj_2.out 366*c4762a1bSJed Brown 367*c4762a1bSJed Brown test: 368*c4762a1bSJed Brown suffix: 6 369*c4762a1bSJed Brown args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_stride 5 -ts_trajectory_solution_only -ts_trajectory_save_stack 0 370*c4762a1bSJed Brown output_file: output/ex20adj_2.out 371*c4762a1bSJed Brown 372*c4762a1bSJed Brown test: 373*c4762a1bSJed Brown suffix: 7 374*c4762a1bSJed Brown args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_stride 5 -ts_trajectory_solution_only 0 -ts_trajectory_save_stack 0 375*c4762a1bSJed Brown output_file: output/ex20adj_2.out 376*c4762a1bSJed Brown 377*c4762a1bSJed Brown test: 378*c4762a1bSJed Brown suffix: 8 379*c4762a1bSJed Brown requires: revolve 380*c4762a1bSJed Brown args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 5 -ts_trajectory_solution_only -ts_trajectory_monitor 381*c4762a1bSJed Brown output_file: output/ex20adj_3.out 382*c4762a1bSJed Brown 383*c4762a1bSJed Brown test: 384*c4762a1bSJed Brown suffix: 9 385*c4762a1bSJed Brown requires: revolve 386*c4762a1bSJed Brown args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 5 -ts_trajectory_solution_only 0 -ts_trajectory_monitor 387*c4762a1bSJed Brown output_file: output/ex20adj_4.out 388*c4762a1bSJed Brown 389*c4762a1bSJed Brown test: 390*c4762a1bSJed Brown requires: revolve 391*c4762a1bSJed Brown suffix: 10 392*c4762a1bSJed Brown args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 5 -ts_trajectory_revolve_online -ts_trajectory_solution_only 393*c4762a1bSJed Brown output_file: output/ex20adj_2.out 394*c4762a1bSJed Brown 395*c4762a1bSJed Brown test: 396*c4762a1bSJed Brown requires: revolve 397*c4762a1bSJed Brown suffix: 11 398*c4762a1bSJed Brown args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 5 -ts_trajectory_revolve_online -ts_trajectory_solution_only 0 399*c4762a1bSJed Brown output_file: output/ex20adj_2.out 400*c4762a1bSJed Brown 401*c4762a1bSJed Brown test: 402*c4762a1bSJed Brown suffix: 12 403*c4762a1bSJed Brown requires: revolve 404*c4762a1bSJed Brown args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 3 -ts_trajectory_max_cps_disk 8 -ts_trajectory_solution_only 405*c4762a1bSJed Brown output_file: output/ex20adj_2.out 406*c4762a1bSJed Brown 407*c4762a1bSJed Brown test: 408*c4762a1bSJed Brown suffix: 13 409*c4762a1bSJed Brown requires: revolve 410*c4762a1bSJed Brown args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 3 -ts_trajectory_max_cps_disk 8 -ts_trajectory_solution_only 0 411*c4762a1bSJed Brown output_file: output/ex20adj_2.out 412*c4762a1bSJed Brown 413*c4762a1bSJed Brown test: 414*c4762a1bSJed Brown suffix: 14 415*c4762a1bSJed Brown requires: revolve 416*c4762a1bSJed Brown args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 3 -ts_trajectory_stride 5 -ts_trajectory_solution_only -ts_trajectory_save_stack 417*c4762a1bSJed Brown output_file: output/ex20adj_2.out 418*c4762a1bSJed Brown 419*c4762a1bSJed Brown test: 420*c4762a1bSJed Brown suffix: 15 421*c4762a1bSJed Brown requires: revolve 422*c4762a1bSJed Brown args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 3 -ts_trajectory_stride 5 -ts_trajectory_solution_only -ts_trajectory_save_stack 0 423*c4762a1bSJed Brown output_file: output/ex20adj_2.out 424*c4762a1bSJed Brown 425*c4762a1bSJed Brown test: 426*c4762a1bSJed Brown suffix: 16 427*c4762a1bSJed Brown requires: revolve 428*c4762a1bSJed Brown args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 3 -ts_trajectory_stride 5 -ts_trajectory_solution_only 0 -ts_trajectory_save_stack 429*c4762a1bSJed Brown output_file: output/ex20adj_2.out 430*c4762a1bSJed Brown 431*c4762a1bSJed Brown test: 432*c4762a1bSJed Brown suffix: 17 433*c4762a1bSJed Brown requires: revolve 434*c4762a1bSJed Brown args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 3 -ts_trajectory_stride 5 -ts_trajectory_solution_only 0 -ts_trajectory_save_stack 0 435*c4762a1bSJed Brown output_file: output/ex20adj_2.out 436*c4762a1bSJed Brown 437*c4762a1bSJed Brown test: 438*c4762a1bSJed Brown suffix: 18 439*c4762a1bSJed Brown requires: revolve 440*c4762a1bSJed Brown args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 3 -ts_trajectory_max_cps_disk 8 -ts_trajectory_stride 5 -ts_trajectory_solution_only -ts_trajectory_save_stack 441*c4762a1bSJed Brown output_file: output/ex20adj_2.out 442*c4762a1bSJed Brown 443*c4762a1bSJed Brown test: 444*c4762a1bSJed Brown suffix: 19 445*c4762a1bSJed Brown requires: revolve 446*c4762a1bSJed Brown args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 3 -ts_trajectory_max_cps_disk 8 -ts_trajectory_stride 5 -ts_trajectory_solution_only 0 -ts_trajectory_save_stack 447*c4762a1bSJed Brown output_file: output/ex20adj_2.out 448*c4762a1bSJed Brown 449*c4762a1bSJed Brown test: 450*c4762a1bSJed Brown suffix: 20 451*c4762a1bSJed Brown requires: revolve 452*c4762a1bSJed Brown args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 3 -ts_trajectory_max_cps_disk 8 -ts_trajectory_solution_only 0 453*c4762a1bSJed Brown output_file: output/ex20adj_2.out 454*c4762a1bSJed Brown 455*c4762a1bSJed Brown test: 456*c4762a1bSJed Brown suffix: 21 457*c4762a1bSJed Brown requires: revolve 458*c4762a1bSJed Brown args: -ts_type cn -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_max_cps_ram 3 -ts_trajectory_max_cps_disk 8 -ts_trajectory_stride 5 -ts_trajectory_solution_only 0 -ts_trajectory_save_stack 0 459*c4762a1bSJed Brown output_file: output/ex20adj_2.out 460*c4762a1bSJed Brown 461*c4762a1bSJed Brown test: 462*c4762a1bSJed Brown suffix: 22 463*c4762a1bSJed Brown args: -ts_type beuler -ts_dt 0.001 -mu 100000 -ts_max_steps 15 -ts_trajectory_type memory -ts_trajectory_solution_only 464*c4762a1bSJed Brown output_file: output/ex20adj_2.out 465*c4762a1bSJed Brown TEST*/ 466