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 REQUIRES configuration of PETSc with option --download-adolc. 7c4762a1bSJed Brown 8c4762a1bSJed Brown For documentation on ADOL-C, see 9c4762a1bSJed Brown $PETSC_ARCH/externalpackages/ADOL-C-2.6.0/ADOL-C/doc/adolc-manual.pdf 10c4762a1bSJed Brown */ 11c4762a1bSJed Brown /* ------------------------------------------------------------------------ 12c4762a1bSJed Brown See ex16adj for a description of the problem being solved. 13c4762a1bSJed Brown ------------------------------------------------------------------------- */ 14c4762a1bSJed Brown 15c4762a1bSJed Brown #include <petscts.h> 16c4762a1bSJed Brown #include <petscmat.h> 17c4762a1bSJed Brown 18c4762a1bSJed Brown #define ADOLC_TAPELESS 19c4762a1bSJed Brown #define NUMBER_DIRECTIONS 3 20c4762a1bSJed Brown #include "adolc-utils/drivers.cxx" 21c4762a1bSJed Brown #include <adolc/adtl.h> 22c4762a1bSJed Brown using namespace adtl; 23c4762a1bSJed Brown 24c4762a1bSJed Brown typedef struct _n_User *User; 25c4762a1bSJed Brown struct _n_User { 26c4762a1bSJed Brown PetscReal mu; 27c4762a1bSJed Brown PetscReal next_output; 28c4762a1bSJed Brown PetscReal tprev; 29c4762a1bSJed Brown 30c4762a1bSJed Brown /* Automatic differentiation support */ 31c4762a1bSJed Brown AdolcCtx *adctx; 32c4762a1bSJed Brown Vec F; 33c4762a1bSJed Brown }; 34c4762a1bSJed Brown 35c4762a1bSJed Brown /* 36c4762a1bSJed Brown Residual evaluation templated, so as to allow for PetscScalar or adouble 37c4762a1bSJed Brown arguments. 38c4762a1bSJed Brown */ 39a8c08197SHong Zhang template <class T> PetscErrorCode EvaluateResidual(const T *x,T mu,T *f) 40c4762a1bSJed Brown { 41c4762a1bSJed Brown PetscFunctionBegin; 42c4762a1bSJed Brown f[0] = x[1]; 43c4762a1bSJed Brown f[1] = mu*(1.-x[0]*x[0])*x[1]-x[0]; 44c4762a1bSJed Brown PetscFunctionReturn(0); 45c4762a1bSJed Brown } 46c4762a1bSJed Brown 47c4762a1bSJed Brown /* 48c4762a1bSJed Brown 'Passive' RHS function, used in residual evaluations during the time integration. 49c4762a1bSJed Brown */ 50c4762a1bSJed Brown static PetscErrorCode RHSFunctionPassive(TS ts,PetscReal t,Vec X,Vec F,void *ctx) 51c4762a1bSJed Brown { 52c4762a1bSJed Brown User user = (User)ctx; 53a8c08197SHong Zhang PetscScalar *f; 54a8c08197SHong Zhang const PetscScalar *x; 55c4762a1bSJed Brown 56c4762a1bSJed Brown PetscFunctionBeginUser; 579566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&x)); 589566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 599566063dSJacob Faibussowitsch PetscCall(EvaluateResidual(x,user->mu,f)); 609566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&x)); 619566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 62c4762a1bSJed Brown PetscFunctionReturn(0); 63c4762a1bSJed Brown } 64c4762a1bSJed Brown 65c4762a1bSJed Brown /* 66c4762a1bSJed Brown Compute the Jacobian w.r.t. x using tapeless mode of ADOL-C. 67c4762a1bSJed Brown */ 68c4762a1bSJed Brown static PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec X,Mat A,Mat B,void *ctx) 69c4762a1bSJed Brown { 70c4762a1bSJed Brown User user = (User)ctx; 71a8c08197SHong Zhang PetscScalar **J; 72a8c08197SHong Zhang const PetscScalar *x; 73c4762a1bSJed Brown adouble f_a[2]; /* 'active' double for dependent variables */ 74c4762a1bSJed Brown adouble x_a[2],mu_a; /* 'active' doubles for independent variables */ 75c4762a1bSJed Brown PetscInt i,j; 76c4762a1bSJed Brown 77c4762a1bSJed Brown PetscFunctionBeginUser; 78c4762a1bSJed Brown /* Set values for independent variables and parameters */ 799566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&x)); 80c4762a1bSJed Brown x_a[0].setValue(x[0]); 81c4762a1bSJed Brown x_a[1].setValue(x[1]); 82c4762a1bSJed Brown mu_a.setValue(user->mu); 839566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&x)); 84c4762a1bSJed Brown 85c4762a1bSJed Brown /* Set seed matrix as 3x3 identity matrix */ 86c4762a1bSJed Brown x_a[0].setADValue(0,1.);x_a[0].setADValue(1,0.);x_a[0].setADValue(2,0.); 87c4762a1bSJed Brown x_a[1].setADValue(0,0.);x_a[1].setADValue(1,1.);x_a[1].setADValue(2,0.); 88c4762a1bSJed Brown mu_a.setADValue(0,0.);mu_a.setADValue(1,0.);mu_a.setADValue(2,1.); 89c4762a1bSJed Brown 90c4762a1bSJed Brown /* Evaluate residual (on active variables) */ 919566063dSJacob Faibussowitsch PetscCall(EvaluateResidual(x_a,mu_a,f_a)); 92c4762a1bSJed Brown 93c4762a1bSJed Brown /* Extract derivatives */ 949566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->adctx->n,&J)); 95c4762a1bSJed Brown J[0] = (PetscScalar*) f_a[0].getADValue(); 96c4762a1bSJed Brown J[1] = (PetscScalar*) f_a[1].getADValue(); 97c4762a1bSJed Brown 98c4762a1bSJed Brown /* Set matrix values */ 99c4762a1bSJed Brown for (i=0; i<user->adctx->m; i++) { 100c4762a1bSJed Brown for (j=0; j<user->adctx->n; j++) { 1019566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&i,1,&j,&J[i][j],INSERT_VALUES)); 102c4762a1bSJed Brown } 103c4762a1bSJed Brown } 1049566063dSJacob Faibussowitsch PetscCall(PetscFree(J)); 1059566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 1069566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 107c4762a1bSJed Brown if (A != B) { 1089566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 1099566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 110c4762a1bSJed Brown } 111c4762a1bSJed Brown PetscFunctionReturn(0); 112c4762a1bSJed Brown } 113c4762a1bSJed Brown 114c4762a1bSJed Brown /* 115c4762a1bSJed Brown Compute the Jacobian w.r.t. mu using tapeless mode of ADOL-C. 116c4762a1bSJed Brown */ 117c4762a1bSJed Brown static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,void *ctx) 118c4762a1bSJed Brown { 119c4762a1bSJed Brown User user = (User)ctx; 120c4762a1bSJed Brown PetscScalar **J; 121c4762a1bSJed Brown PetscScalar *x; 122c4762a1bSJed Brown adouble f_a[2]; /* 'active' double for dependent variables */ 123c4762a1bSJed Brown adouble x_a[2],mu_a; /* 'active' doubles for independent variables */ 124c4762a1bSJed Brown PetscInt i,j = 0; 125c4762a1bSJed Brown 126c4762a1bSJed Brown PetscFunctionBeginUser; 127c4762a1bSJed Brown 128c4762a1bSJed Brown /* Set values for independent variables and parameters */ 1299566063dSJacob Faibussowitsch PetscCall(VecGetArray(X,&x)); 130c4762a1bSJed Brown x_a[0].setValue(x[0]); 131c4762a1bSJed Brown x_a[1].setValue(x[1]); 132c4762a1bSJed Brown mu_a.setValue(user->mu); 1339566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(X,&x)); 134c4762a1bSJed Brown 135c4762a1bSJed Brown /* Set seed matrix as 3x3 identity matrix */ 136c4762a1bSJed Brown x_a[0].setADValue(0,1.);x_a[0].setADValue(1,0.);x_a[0].setADValue(2,0.); 137c4762a1bSJed Brown x_a[1].setADValue(0,0.);x_a[1].setADValue(1,1.);x_a[1].setADValue(2,0.); 138c4762a1bSJed Brown mu_a.setADValue(0,0.);mu_a.setADValue(1,0.);mu_a.setADValue(2,1.); 139c4762a1bSJed Brown 140c4762a1bSJed Brown /* Evaluate residual (on active variables) */ 1419566063dSJacob Faibussowitsch PetscCall(EvaluateResidual(x_a,mu_a,f_a)); 142c4762a1bSJed Brown 143c4762a1bSJed Brown /* Extract derivatives */ 1449566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2,&J)); 145c4762a1bSJed Brown J[0] = (PetscScalar*) f_a[0].getADValue(); 146c4762a1bSJed Brown J[1] = (PetscScalar*) f_a[1].getADValue(); 147c4762a1bSJed Brown 148c4762a1bSJed Brown /* Set matrix values */ 149c4762a1bSJed Brown for (i=0; i<user->adctx->m; i++) { 1509566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,&i,1,&j,&J[i][user->adctx->n],INSERT_VALUES)); 151c4762a1bSJed Brown } 1529566063dSJacob Faibussowitsch PetscCall(PetscFree(J)); 1539566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 1549566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 155c4762a1bSJed Brown PetscFunctionReturn(0); 156c4762a1bSJed Brown } 157c4762a1bSJed Brown 158c4762a1bSJed Brown /* 159c4762a1bSJed Brown Monitor timesteps and use interpolation to output at integer multiples of 0.1 160c4762a1bSJed Brown */ 161c4762a1bSJed Brown static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx) 162c4762a1bSJed Brown { 163c4762a1bSJed Brown const PetscScalar *x; 164c4762a1bSJed Brown PetscReal tfinal, dt, tprev; 165c4762a1bSJed Brown User user = (User)ctx; 166c4762a1bSJed Brown 167c4762a1bSJed Brown PetscFunctionBeginUser; 1689566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts,&dt)); 1699566063dSJacob Faibussowitsch PetscCall(TSGetMaxTime(ts,&tfinal)); 1709566063dSJacob Faibussowitsch PetscCall(TSGetPrevTime(ts,&tprev)); 1719566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&x)); 17263a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"[%.1f] %" PetscInt_FMT " 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]))); 1739566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"t %.6f (tprev = %.6f) \n",(double)t,(double)tprev)); 1749566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&x)); 175c4762a1bSJed Brown PetscFunctionReturn(0); 176c4762a1bSJed Brown } 177c4762a1bSJed Brown 178c4762a1bSJed Brown int main(int argc,char **argv) 179c4762a1bSJed Brown { 180c4762a1bSJed Brown TS ts; /* nonlinear solver */ 181c4762a1bSJed Brown Vec x; /* solution, residual vectors */ 182c4762a1bSJed Brown Mat A; /* Jacobian matrix */ 183c4762a1bSJed Brown Mat Jacp; /* JacobianP matrix */ 184c4762a1bSJed Brown PetscInt steps; 185c4762a1bSJed Brown PetscReal ftime = 0.5; 186c4762a1bSJed Brown PetscBool monitor = PETSC_FALSE; 187c4762a1bSJed Brown PetscScalar *x_ptr; 188c4762a1bSJed Brown PetscMPIInt size; 189c4762a1bSJed Brown struct _n_User user; 190c4762a1bSJed Brown AdolcCtx *adctx; 191c4762a1bSJed Brown Vec lambda[2],mu[2]; 192c4762a1bSJed Brown 193c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 194c4762a1bSJed Brown Initialize program 195c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 196*327415f7SBarry Smith PetscFunctionBeginUser; 1979566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,NULL,help)); 1989566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 19908401ef6SPierre Jolivet PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 200c4762a1bSJed Brown 201c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 202c4762a1bSJed Brown Set runtime options and create AdolcCtx 203c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2049566063dSJacob Faibussowitsch PetscCall(PetscNew(&adctx)); 205c4762a1bSJed Brown user.mu = 1; 206c4762a1bSJed Brown user.next_output = 0.0; 207c4762a1bSJed Brown adctx->m = 2;adctx->n = 2;adctx->p = 2; 208c4762a1bSJed Brown user.adctx = adctx; 209c4762a1bSJed Brown adtl::setNumDir(adctx->n+1); /* #indep. variables, plus parameters */ 210c4762a1bSJed Brown 2119566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,NULL)); 2129566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL)); 213c4762a1bSJed Brown 214c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 215c4762a1bSJed Brown Create necessary matrix and vectors, solve same ODE on every process 216c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2179566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); 2189566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,2,2)); 2199566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 2209566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 2219566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A,&x,NULL)); 222c4762a1bSJed Brown 2239566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&Jacp)); 2249566063dSJacob Faibussowitsch PetscCall(MatSetSizes(Jacp,PETSC_DECIDE,PETSC_DECIDE,2,1)); 2259566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(Jacp)); 2269566063dSJacob Faibussowitsch PetscCall(MatSetUp(Jacp)); 227c4762a1bSJed Brown 228c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 229c4762a1bSJed Brown Create timestepping solver context 230c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2319566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD,&ts)); 2329566063dSJacob Faibussowitsch PetscCall(TSSetType(ts,TSRK)); 2339566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts,NULL,RHSFunctionPassive,&user)); 234c4762a1bSJed Brown 235c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 236c4762a1bSJed Brown Set initial conditions 237c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2389566063dSJacob Faibussowitsch PetscCall(VecGetArray(x,&x_ptr)); 239c4762a1bSJed Brown x_ptr[0] = 2; x_ptr[1] = 0.66666654321; 2409566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x,&x_ptr)); 241c4762a1bSJed Brown 242c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 243c4762a1bSJed Brown Set RHS Jacobian for the adjoint integration 244c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2459566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(ts,A,A,RHSJacobian,&user)); 2469566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts,ftime)); 2479566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP)); 248c4762a1bSJed Brown if (monitor) { 2499566063dSJacob Faibussowitsch PetscCall(TSMonitorSet(ts,Monitor,&user,NULL)); 250c4762a1bSJed Brown } 2519566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts,.001)); 252c4762a1bSJed Brown 253c4762a1bSJed Brown /* 254c4762a1bSJed Brown Have the TS save its trajectory so that TSAdjointSolve() may be used 255c4762a1bSJed Brown */ 2569566063dSJacob Faibussowitsch PetscCall(TSSetSaveTrajectory(ts)); 257c4762a1bSJed Brown 258c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 259c4762a1bSJed Brown Set runtime options 260c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2619566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 262c4762a1bSJed Brown 263c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 264c4762a1bSJed Brown Solve nonlinear system 265c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2669566063dSJacob Faibussowitsch PetscCall(TSSolve(ts,x)); 2679566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts,&ftime)); 2689566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts,&steps)); 26963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %" PetscInt_FMT ", ftime %g\n",(double)user.mu,steps,(double)ftime)); 2709566063dSJacob Faibussowitsch PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD)); 271c4762a1bSJed Brown 272c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 273c4762a1bSJed Brown Start the Adjoint model 274c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2759566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A,&lambda[0],NULL)); 2769566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A,&lambda[1],NULL)); 277c4762a1bSJed Brown /* Reset initial conditions for the adjoint integration */ 2789566063dSJacob Faibussowitsch PetscCall(VecGetArray(lambda[0],&x_ptr)); 279c4762a1bSJed Brown x_ptr[0] = 1.0; x_ptr[1] = 0.0; 2809566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(lambda[0],&x_ptr)); 2819566063dSJacob Faibussowitsch PetscCall(VecGetArray(lambda[1],&x_ptr)); 282c4762a1bSJed Brown x_ptr[0] = 0.0; x_ptr[1] = 1.0; 2839566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(lambda[1],&x_ptr)); 284c4762a1bSJed Brown 2859566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(Jacp,&mu[0],NULL)); 2869566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(Jacp,&mu[1],NULL)); 2879566063dSJacob Faibussowitsch PetscCall(VecGetArray(mu[0],&x_ptr)); 288c4762a1bSJed Brown x_ptr[0] = 0.0; 2899566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(mu[0],&x_ptr)); 2909566063dSJacob Faibussowitsch PetscCall(VecGetArray(mu[1],&x_ptr)); 291c4762a1bSJed Brown x_ptr[0] = 0.0; 2929566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(mu[1],&x_ptr)); 2939566063dSJacob Faibussowitsch PetscCall(TSSetCostGradients(ts,2,lambda,mu)); 294c4762a1bSJed Brown 295c4762a1bSJed Brown /* Set RHS JacobianP */ 2969566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobianP(ts,Jacp,RHSJacobianP,&user)); 297c4762a1bSJed Brown 2989566063dSJacob Faibussowitsch PetscCall(TSAdjointSolve(ts)); 299c4762a1bSJed Brown 3009566063dSJacob Faibussowitsch PetscCall(VecView(lambda[0],PETSC_VIEWER_STDOUT_WORLD)); 3019566063dSJacob Faibussowitsch PetscCall(VecView(lambda[1],PETSC_VIEWER_STDOUT_WORLD)); 3029566063dSJacob Faibussowitsch PetscCall(VecView(mu[0],PETSC_VIEWER_STDOUT_WORLD)); 3039566063dSJacob Faibussowitsch PetscCall(VecView(mu[1],PETSC_VIEWER_STDOUT_WORLD)); 304c4762a1bSJed Brown 305c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 306c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 307c4762a1bSJed Brown are no longer needed. 308c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 3099566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 3109566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Jacp)); 3119566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 3129566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lambda[0])); 3139566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lambda[1])); 3149566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mu[0])); 3159566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mu[1])); 3169566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 3179566063dSJacob Faibussowitsch PetscCall(PetscFree(adctx)); 3189566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 319b122ec5aSJacob Faibussowitsch return 0; 320c4762a1bSJed Brown } 321c4762a1bSJed Brown 322c4762a1bSJed Brown /*TEST 323c4762a1bSJed Brown 324c4762a1bSJed Brown build: 325c4762a1bSJed Brown requires: double !complex adolc 326c4762a1bSJed Brown 327c4762a1bSJed Brown test: 328c4762a1bSJed Brown suffix: 1 329c4762a1bSJed Brown args: -ts_max_steps 10 -ts_monitor -ts_adjoint_monitor 330c4762a1bSJed Brown output_file: output/ex16adj_tl_1.out 331c4762a1bSJed Brown 332c4762a1bSJed Brown test: 333c4762a1bSJed Brown suffix: 2 334c4762a1bSJed Brown args: -ts_max_steps 10 -ts_monitor -ts_adjoint_monitor -mu 5 335c4762a1bSJed Brown output_file: output/ex16adj_tl_2.out 336c4762a1bSJed Brown 337c4762a1bSJed Brown test: 338c4762a1bSJed Brown suffix: 3 339c4762a1bSJed Brown args: -ts_max_steps 10 -monitor 340c4762a1bSJed Brown output_file: output/ex16adj_tl_3.out 341c4762a1bSJed Brown 342c4762a1bSJed Brown test: 343c4762a1bSJed Brown suffix: 4 344c4762a1bSJed Brown args: -ts_max_steps 10 -monitor -mu 5 345c4762a1bSJed Brown output_file: output/ex16adj_tl_4.out 346c4762a1bSJed Brown 347c4762a1bSJed Brown TEST*/ 348