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