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