1*c4762a1bSJed Brown static char help[] = "Demonstrates automatic Jacobian generation using ADOL-C for an adjoint sensitivity analysis of the van der Pol equation.\n\ 2*c4762a1bSJed Brown Input parameters include:\n\ 3*c4762a1bSJed Brown -mu : stiffness parameter\n\n"; 4*c4762a1bSJed Brown 5*c4762a1bSJed Brown /* 6*c4762a1bSJed Brown Concepts: TS^time-dependent nonlinear problems 7*c4762a1bSJed Brown Concepts: TS^van der Pol equation 8*c4762a1bSJed Brown Concepts: TS^adjoint sensitivity analysis 9*c4762a1bSJed Brown Concepts: Automatic differentation using ADOL-C 10*c4762a1bSJed Brown Concepts: Automatic differentation w.r.t. a parameter using ADOL-C 11*c4762a1bSJed Brown Processors: 1 12*c4762a1bSJed Brown */ 13*c4762a1bSJed Brown /* 14*c4762a1bSJed Brown REQUIRES configuration of PETSc with option --download-adolc. 15*c4762a1bSJed Brown 16*c4762a1bSJed Brown For documentation on ADOL-C, see 17*c4762a1bSJed Brown $PETSC_ARCH/externalpackages/ADOL-C-2.6.0/ADOL-C/doc/adolc-manual.pdf 18*c4762a1bSJed Brown */ 19*c4762a1bSJed Brown /* ------------------------------------------------------------------------ 20*c4762a1bSJed Brown See ex16adj for a description of the problem being solved. 21*c4762a1bSJed Brown ------------------------------------------------------------------------- */ 22*c4762a1bSJed Brown 23*c4762a1bSJed Brown #include <petscts.h> 24*c4762a1bSJed Brown #include <petscmat.h> 25*c4762a1bSJed Brown #include "adolc-utils/drivers.cxx" 26*c4762a1bSJed Brown #include <adolc/adolc.h> 27*c4762a1bSJed Brown 28*c4762a1bSJed Brown typedef struct _n_User *User; 29*c4762a1bSJed Brown struct _n_User { 30*c4762a1bSJed Brown PetscReal mu; 31*c4762a1bSJed Brown PetscReal next_output; 32*c4762a1bSJed Brown PetscReal tprev; 33*c4762a1bSJed Brown 34*c4762a1bSJed Brown /* Automatic differentiation support */ 35*c4762a1bSJed Brown AdolcCtx *adctx; 36*c4762a1bSJed Brown }; 37*c4762a1bSJed Brown 38*c4762a1bSJed Brown /* 39*c4762a1bSJed Brown 'Passive' RHS function, used in residual evaluations during the time integration. 40*c4762a1bSJed Brown */ 41*c4762a1bSJed Brown static PetscErrorCode RHSFunctionPassive(TS ts,PetscReal t,Vec X,Vec F,void *ctx) 42*c4762a1bSJed Brown { 43*c4762a1bSJed Brown PetscErrorCode ierr; 44*c4762a1bSJed Brown User user = (User)ctx; 45*c4762a1bSJed Brown PetscScalar *f; 46*c4762a1bSJed Brown const PetscScalar *x; 47*c4762a1bSJed Brown 48*c4762a1bSJed Brown PetscFunctionBeginUser; 49*c4762a1bSJed Brown ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 50*c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 51*c4762a1bSJed Brown f[0] = x[1]; 52*c4762a1bSJed Brown f[1] = user->mu*(1.-x[0]*x[0])*x[1]-x[0]; 53*c4762a1bSJed Brown ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 54*c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 55*c4762a1bSJed Brown PetscFunctionReturn(0); 56*c4762a1bSJed Brown } 57*c4762a1bSJed Brown 58*c4762a1bSJed Brown /* 59*c4762a1bSJed Brown Trace RHS to mark on tape 1 the dependence of f upon x. This tape is used in generating the 60*c4762a1bSJed Brown Jacobian transform. 61*c4762a1bSJed Brown */ 62*c4762a1bSJed Brown static PetscErrorCode RHSFunctionActive(TS ts,PetscReal t,Vec X,Vec F,void *ctx) 63*c4762a1bSJed Brown { 64*c4762a1bSJed Brown PetscErrorCode ierr; 65*c4762a1bSJed Brown User user = (User)ctx; 66*c4762a1bSJed Brown PetscScalar *f; 67*c4762a1bSJed Brown const PetscScalar *x; 68*c4762a1bSJed Brown 69*c4762a1bSJed Brown adouble f_a[2]; /* 'active' double for dependent variables */ 70*c4762a1bSJed Brown adouble x_a[2]; /* 'active' double for independent variables */ 71*c4762a1bSJed Brown 72*c4762a1bSJed Brown PetscFunctionBeginUser; 73*c4762a1bSJed Brown ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 74*c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 75*c4762a1bSJed Brown 76*c4762a1bSJed Brown /* Start of active section */ 77*c4762a1bSJed Brown trace_on(1); 78*c4762a1bSJed Brown x_a[0] <<= x[0];x_a[1] <<= x[1]; /* Mark independence */ 79*c4762a1bSJed Brown f_a[0] = x_a[1]; 80*c4762a1bSJed Brown f_a[1] = user->mu*(1.-x_a[0]*x_a[0])*x_a[1]-x_a[0]; 81*c4762a1bSJed Brown f_a[0] >>= f[0];f_a[1] >>= f[1]; /* Mark dependence */ 82*c4762a1bSJed Brown trace_off(); 83*c4762a1bSJed Brown /* End of active section */ 84*c4762a1bSJed Brown 85*c4762a1bSJed Brown ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 86*c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 87*c4762a1bSJed Brown PetscFunctionReturn(0); 88*c4762a1bSJed Brown } 89*c4762a1bSJed Brown 90*c4762a1bSJed Brown /* 91*c4762a1bSJed Brown Trace RHS again to mark on tape 2 the dependence of f upon the parameter mu. This tape is used in 92*c4762a1bSJed Brown generating JacobianP. 93*c4762a1bSJed Brown */ 94*c4762a1bSJed Brown static PetscErrorCode RHSFunctionActiveP(TS ts,PetscReal t,Vec X,Vec F,void *ctx) 95*c4762a1bSJed Brown { 96*c4762a1bSJed Brown PetscErrorCode ierr; 97*c4762a1bSJed Brown User user = (User)ctx; 98*c4762a1bSJed Brown PetscScalar *f; 99*c4762a1bSJed Brown const PetscScalar *x; 100*c4762a1bSJed Brown 101*c4762a1bSJed Brown adouble f_a[2]; /* 'active' double for dependent variables */ 102*c4762a1bSJed Brown adouble x_a[2],mu_a; /* 'active' double for independent variables */ 103*c4762a1bSJed Brown 104*c4762a1bSJed Brown PetscFunctionBeginUser; 105*c4762a1bSJed Brown ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 106*c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 107*c4762a1bSJed Brown 108*c4762a1bSJed Brown /* Start of active section */ 109*c4762a1bSJed Brown trace_on(3); 110*c4762a1bSJed Brown x_a[0] <<= x[0];x_a[1] <<= x[1];mu_a <<= user->mu; /* Mark independence */ 111*c4762a1bSJed Brown f_a[0] = x_a[1]; 112*c4762a1bSJed Brown f_a[1] = mu_a*(1.-x_a[0]*x_a[0])*x_a[1]-x_a[0]; 113*c4762a1bSJed Brown f_a[0] >>= f[0];f_a[1] >>= f[1]; /* Mark dependence */ 114*c4762a1bSJed Brown trace_off(); 115*c4762a1bSJed Brown /* End of active section */ 116*c4762a1bSJed Brown 117*c4762a1bSJed Brown ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 118*c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 119*c4762a1bSJed Brown PetscFunctionReturn(0); 120*c4762a1bSJed Brown } 121*c4762a1bSJed Brown 122*c4762a1bSJed Brown /* 123*c4762a1bSJed Brown Compute the Jacobian w.r.t. x using PETSc-ADOL-C driver for explicit TS. 124*c4762a1bSJed Brown */ 125*c4762a1bSJed Brown static PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec X,Mat A,Mat B,void *ctx) 126*c4762a1bSJed Brown { 127*c4762a1bSJed Brown PetscErrorCode ierr; 128*c4762a1bSJed Brown User user = (User)ctx; 129*c4762a1bSJed Brown PetscScalar *x; 130*c4762a1bSJed Brown 131*c4762a1bSJed Brown PetscFunctionBeginUser; 132*c4762a1bSJed Brown ierr = VecGetArray(X,&x);CHKERRQ(ierr); 133*c4762a1bSJed Brown ierr = PetscAdolcComputeRHSJacobian(1,A,x,user->adctx);CHKERRQ(ierr); 134*c4762a1bSJed Brown ierr = VecRestoreArray(X,&x);CHKERRQ(ierr); 135*c4762a1bSJed Brown PetscFunctionReturn(0); 136*c4762a1bSJed Brown } 137*c4762a1bSJed Brown 138*c4762a1bSJed Brown /* 139*c4762a1bSJed Brown Compute the Jacobian w.r.t. mu using PETSc-ADOL-C driver for explicit TS. 140*c4762a1bSJed Brown */ 141*c4762a1bSJed Brown static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,void *ctx) 142*c4762a1bSJed Brown { 143*c4762a1bSJed Brown PetscErrorCode ierr; 144*c4762a1bSJed Brown User user = (User)ctx; 145*c4762a1bSJed Brown PetscScalar *x; 146*c4762a1bSJed Brown 147*c4762a1bSJed Brown PetscFunctionBeginUser; 148*c4762a1bSJed Brown ierr = VecGetArray(X,&x);CHKERRQ(ierr); 149*c4762a1bSJed Brown ierr = PetscAdolcComputeRHSJacobianP(3,A,x,&user->mu,user->adctx);CHKERRQ(ierr); 150*c4762a1bSJed Brown ierr = VecRestoreArray(X,&x);CHKERRQ(ierr); 151*c4762a1bSJed Brown PetscFunctionReturn(0); 152*c4762a1bSJed Brown } 153*c4762a1bSJed Brown 154*c4762a1bSJed Brown /* 155*c4762a1bSJed Brown Monitor timesteps and use interpolation to output at integer multiples of 0.1 156*c4762a1bSJed Brown */ 157*c4762a1bSJed Brown static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx) 158*c4762a1bSJed Brown { 159*c4762a1bSJed Brown PetscErrorCode ierr; 160*c4762a1bSJed Brown const PetscScalar *x; 161*c4762a1bSJed Brown PetscReal tfinal, dt, tprev; 162*c4762a1bSJed Brown User user = (User)ctx; 163*c4762a1bSJed Brown 164*c4762a1bSJed Brown PetscFunctionBeginUser; 165*c4762a1bSJed Brown ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr); 166*c4762a1bSJed Brown ierr = TSGetMaxTime(ts,&tfinal);CHKERRQ(ierr); 167*c4762a1bSJed Brown ierr = TSGetPrevTime(ts,&tprev);CHKERRQ(ierr); 168*c4762a1bSJed Brown ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 169*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"[%.1f] %D TS %.6f (dt = %.6f) X % 12.6e % 12.6e\n",(double)user->next_output,step,(double)t,(double)dt,(double)PetscRealPart(x[0]),(double)PetscRealPart(x[1]));CHKERRQ(ierr); 170*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"t %.6f (tprev = %.6f) \n",(double)t,(double)tprev);CHKERRQ(ierr); 171*c4762a1bSJed Brown ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 172*c4762a1bSJed Brown PetscFunctionReturn(0); 173*c4762a1bSJed Brown } 174*c4762a1bSJed Brown 175*c4762a1bSJed Brown int main(int argc,char **argv) 176*c4762a1bSJed Brown { 177*c4762a1bSJed Brown TS ts; /* nonlinear solver */ 178*c4762a1bSJed Brown Vec x; /* solution, residual vectors */ 179*c4762a1bSJed Brown Mat A; /* Jacobian matrix */ 180*c4762a1bSJed Brown Mat Jacp; /* JacobianP matrix */ 181*c4762a1bSJed Brown PetscInt steps; 182*c4762a1bSJed Brown PetscReal ftime = 0.5; 183*c4762a1bSJed Brown PetscBool monitor = PETSC_FALSE; 184*c4762a1bSJed Brown PetscScalar *x_ptr; 185*c4762a1bSJed Brown PetscMPIInt size; 186*c4762a1bSJed Brown struct _n_User user; 187*c4762a1bSJed Brown AdolcCtx *adctx; 188*c4762a1bSJed Brown PetscErrorCode ierr; 189*c4762a1bSJed Brown Vec lambda[2],mu[2],r; 190*c4762a1bSJed Brown 191*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 192*c4762a1bSJed Brown Initialize program 193*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 194*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr; 195*c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 196*c4762a1bSJed Brown if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_MPI_WRONG_SIZE,"This is a uniprocessor example only!"); 197*c4762a1bSJed Brown 198*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 199*c4762a1bSJed Brown Set runtime options and create AdolcCtx 200*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 201*c4762a1bSJed Brown ierr = PetscNew(&adctx);CHKERRQ(ierr); 202*c4762a1bSJed Brown user.mu = 1; 203*c4762a1bSJed Brown user.next_output = 0.0; 204*c4762a1bSJed Brown adctx->m = 2;adctx->n = 2;adctx->p = 2;adctx->num_params = 1; 205*c4762a1bSJed Brown user.adctx = adctx; 206*c4762a1bSJed Brown 207*c4762a1bSJed Brown ierr = PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,NULL);CHKERRQ(ierr); 208*c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL);CHKERRQ(ierr); 209*c4762a1bSJed Brown 210*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 211*c4762a1bSJed Brown Create necessary matrix and vectors, solve same ODE on every process 212*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 213*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 214*c4762a1bSJed Brown ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,2,2);CHKERRQ(ierr); 215*c4762a1bSJed Brown ierr = MatSetFromOptions(A);CHKERRQ(ierr); 216*c4762a1bSJed Brown ierr = MatSetUp(A);CHKERRQ(ierr); 217*c4762a1bSJed Brown ierr = MatCreateVecs(A,&x,NULL);CHKERRQ(ierr); 218*c4762a1bSJed Brown 219*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&Jacp);CHKERRQ(ierr); 220*c4762a1bSJed Brown ierr = MatSetSizes(Jacp,PETSC_DECIDE,PETSC_DECIDE,2,1);CHKERRQ(ierr); 221*c4762a1bSJed Brown ierr = MatSetFromOptions(Jacp);CHKERRQ(ierr); 222*c4762a1bSJed Brown ierr = MatSetUp(Jacp);CHKERRQ(ierr); 223*c4762a1bSJed Brown 224*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 225*c4762a1bSJed Brown Create timestepping solver context 226*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 227*c4762a1bSJed Brown ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr); 228*c4762a1bSJed Brown ierr = TSSetType(ts,TSRK);CHKERRQ(ierr); 229*c4762a1bSJed Brown ierr = TSSetRHSFunction(ts,NULL,RHSFunctionPassive,&user);CHKERRQ(ierr); 230*c4762a1bSJed Brown 231*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 232*c4762a1bSJed Brown Set initial conditions 233*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 234*c4762a1bSJed Brown ierr = VecGetArray(x,&x_ptr);CHKERRQ(ierr); 235*c4762a1bSJed Brown x_ptr[0] = 2; x_ptr[1] = 0.66666654321; 236*c4762a1bSJed Brown ierr = VecRestoreArray(x,&x_ptr);CHKERRQ(ierr); 237*c4762a1bSJed Brown 238*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 239*c4762a1bSJed Brown Trace just once on each tape and put zeros on Jacobian diagonal 240*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 241*c4762a1bSJed Brown ierr = VecDuplicate(x,&r);CHKERRQ(ierr); 242*c4762a1bSJed Brown ierr = RHSFunctionActive(ts,0.,x,r,&user);CHKERRQ(ierr); 243*c4762a1bSJed Brown ierr = RHSFunctionActiveP(ts,0.,x,r,&user);CHKERRQ(ierr); 244*c4762a1bSJed Brown ierr = VecSet(r,0);CHKERRQ(ierr); 245*c4762a1bSJed Brown ierr = MatDiagonalSet(A,r,INSERT_VALUES);CHKERRQ(ierr); 246*c4762a1bSJed Brown ierr = VecDestroy(&r);CHKERRQ(ierr); 247*c4762a1bSJed Brown 248*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 249*c4762a1bSJed Brown Set RHS Jacobian for the adjoint integration 250*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 251*c4762a1bSJed Brown ierr = TSSetRHSJacobian(ts,A,A,RHSJacobian,&user);CHKERRQ(ierr); 252*c4762a1bSJed Brown ierr = TSSetMaxTime(ts,ftime);CHKERRQ(ierr); 253*c4762a1bSJed Brown ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr); 254*c4762a1bSJed Brown if (monitor) { 255*c4762a1bSJed Brown ierr = TSMonitorSet(ts,Monitor,&user,NULL);CHKERRQ(ierr); 256*c4762a1bSJed Brown } 257*c4762a1bSJed Brown ierr = TSSetTimeStep(ts,.001);CHKERRQ(ierr); 258*c4762a1bSJed Brown 259*c4762a1bSJed Brown /* 260*c4762a1bSJed Brown Have the TS save its trajectory so that TSAdjointSolve() may be used 261*c4762a1bSJed Brown */ 262*c4762a1bSJed Brown ierr = TSSetSaveTrajectory(ts);CHKERRQ(ierr); 263*c4762a1bSJed Brown 264*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 265*c4762a1bSJed Brown Set runtime options 266*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 267*c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 268*c4762a1bSJed Brown 269*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 270*c4762a1bSJed Brown Solve nonlinear system 271*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 272*c4762a1bSJed Brown ierr = TSSolve(ts,x);CHKERRQ(ierr); 273*c4762a1bSJed Brown ierr = TSGetSolveTime(ts,&ftime);CHKERRQ(ierr); 274*c4762a1bSJed Brown ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr); 275*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %D, ftime %g\n",(double)user.mu,steps,(double)ftime);CHKERRQ(ierr); 276*c4762a1bSJed Brown ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 277*c4762a1bSJed Brown 278*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 279*c4762a1bSJed Brown Start the Adjoint model 280*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 281*c4762a1bSJed Brown ierr = MatCreateVecs(A,&lambda[0],NULL);CHKERRQ(ierr); 282*c4762a1bSJed Brown ierr = MatCreateVecs(A,&lambda[1],NULL);CHKERRQ(ierr); 283*c4762a1bSJed Brown /* Reset initial conditions for the adjoint integration */ 284*c4762a1bSJed Brown ierr = VecGetArray(lambda[0],&x_ptr);CHKERRQ(ierr); 285*c4762a1bSJed Brown x_ptr[0] = 1.0; x_ptr[1] = 0.0; 286*c4762a1bSJed Brown ierr = VecRestoreArray(lambda[0],&x_ptr);CHKERRQ(ierr); 287*c4762a1bSJed Brown ierr = VecGetArray(lambda[1],&x_ptr);CHKERRQ(ierr); 288*c4762a1bSJed Brown x_ptr[0] = 0.0; x_ptr[1] = 1.0; 289*c4762a1bSJed Brown ierr = VecRestoreArray(lambda[1],&x_ptr);CHKERRQ(ierr); 290*c4762a1bSJed Brown 291*c4762a1bSJed Brown ierr = MatCreateVecs(Jacp,&mu[0],NULL);CHKERRQ(ierr); 292*c4762a1bSJed Brown ierr = MatCreateVecs(Jacp,&mu[1],NULL);CHKERRQ(ierr); 293*c4762a1bSJed Brown ierr = VecGetArray(mu[0],&x_ptr);CHKERRQ(ierr); 294*c4762a1bSJed Brown x_ptr[0] = 0.0; 295*c4762a1bSJed Brown ierr = VecRestoreArray(mu[0],&x_ptr);CHKERRQ(ierr); 296*c4762a1bSJed Brown ierr = VecGetArray(mu[1],&x_ptr);CHKERRQ(ierr); 297*c4762a1bSJed Brown x_ptr[0] = 0.0; 298*c4762a1bSJed Brown ierr = VecRestoreArray(mu[1],&x_ptr);CHKERRQ(ierr); 299*c4762a1bSJed Brown ierr = TSSetCostGradients(ts,2,lambda,mu);CHKERRQ(ierr); 300*c4762a1bSJed Brown 301*c4762a1bSJed Brown 302*c4762a1bSJed Brown /* Set RHS JacobianP */ 303*c4762a1bSJed Brown ierr = TSSetRHSJacobianP(ts,Jacp,RHSJacobianP,&user);CHKERRQ(ierr); 304*c4762a1bSJed Brown 305*c4762a1bSJed Brown ierr = TSAdjointSolve(ts);CHKERRQ(ierr); 306*c4762a1bSJed Brown 307*c4762a1bSJed Brown ierr = VecView(lambda[0],PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 308*c4762a1bSJed Brown ierr = VecView(lambda[1],PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 309*c4762a1bSJed Brown ierr = VecView(mu[0],PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 310*c4762a1bSJed Brown ierr = VecView(mu[1],PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 311*c4762a1bSJed Brown 312*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 313*c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 314*c4762a1bSJed Brown are no longer needed. 315*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 316*c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 317*c4762a1bSJed Brown ierr = MatDestroy(&Jacp);CHKERRQ(ierr); 318*c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); 319*c4762a1bSJed Brown ierr = VecDestroy(&lambda[0]);CHKERRQ(ierr); 320*c4762a1bSJed Brown ierr = VecDestroy(&lambda[1]);CHKERRQ(ierr); 321*c4762a1bSJed Brown ierr = VecDestroy(&mu[0]);CHKERRQ(ierr); 322*c4762a1bSJed Brown ierr = VecDestroy(&mu[1]);CHKERRQ(ierr); 323*c4762a1bSJed Brown ierr = TSDestroy(&ts);CHKERRQ(ierr); 324*c4762a1bSJed Brown ierr = PetscFree(adctx);CHKERRQ(ierr); 325*c4762a1bSJed Brown ierr = PetscFinalize(); 326*c4762a1bSJed Brown return ierr; 327*c4762a1bSJed Brown } 328*c4762a1bSJed Brown 329*c4762a1bSJed Brown /*TEST 330*c4762a1bSJed Brown 331*c4762a1bSJed Brown build: 332*c4762a1bSJed Brown requires: double !complex adolc 333*c4762a1bSJed Brown 334*c4762a1bSJed Brown test: 335*c4762a1bSJed Brown suffix: 1 336*c4762a1bSJed Brown args: -ts_max_steps 10 -ts_monitor -ts_adjoint_monitor 337*c4762a1bSJed Brown output_file: output/ex16adj_1.out 338*c4762a1bSJed Brown 339*c4762a1bSJed Brown test: 340*c4762a1bSJed Brown suffix: 2 341*c4762a1bSJed Brown args: -ts_max_steps 10 -ts_monitor -ts_adjoint_monitor -mu 5 342*c4762a1bSJed Brown output_file: output/ex16adj_2.out 343*c4762a1bSJed Brown 344*c4762a1bSJed Brown test: 345*c4762a1bSJed Brown suffix: 3 346*c4762a1bSJed Brown args: -ts_max_steps 10 -monitor 347*c4762a1bSJed Brown output_file: output/ex16adj_3.out 348*c4762a1bSJed Brown 349*c4762a1bSJed Brown test: 350*c4762a1bSJed Brown suffix: 4 351*c4762a1bSJed Brown args: -ts_max_steps 10 -monitor -mu 5 352*c4762a1bSJed Brown output_file: output/ex16adj_4.out 353*c4762a1bSJed Brown 354*c4762a1bSJed Brown TEST*/ 355