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 */ 399371c9d4SSatish Balay template <class T> 409371c9d4SSatish Balay PetscErrorCode EvaluateResidual(const T *x, T mu, T *f) { 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 */ 509371c9d4SSatish Balay static PetscErrorCode RHSFunctionPassive(TS ts, PetscReal t, Vec X, Vec F, void *ctx) { 51c4762a1bSJed Brown User user = (User)ctx; 52a8c08197SHong Zhang PetscScalar *f; 53a8c08197SHong Zhang const PetscScalar *x; 54c4762a1bSJed Brown 55c4762a1bSJed Brown PetscFunctionBeginUser; 569566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 579566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 589566063dSJacob Faibussowitsch PetscCall(EvaluateResidual(x, user->mu, f)); 599566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 609566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 61c4762a1bSJed Brown PetscFunctionReturn(0); 62c4762a1bSJed Brown } 63c4762a1bSJed Brown 64c4762a1bSJed Brown /* 65c4762a1bSJed Brown Compute the Jacobian w.r.t. x using tapeless mode of ADOL-C. 66c4762a1bSJed Brown */ 679371c9d4SSatish Balay static PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec X, Mat A, Mat B, void *ctx) { 68c4762a1bSJed Brown User user = (User)ctx; 69a8c08197SHong Zhang PetscScalar **J; 70a8c08197SHong Zhang const PetscScalar *x; 71c4762a1bSJed Brown adouble f_a[2]; /* 'active' double for dependent variables */ 72c4762a1bSJed Brown adouble x_a[2], mu_a; /* 'active' doubles for independent variables */ 73c4762a1bSJed Brown PetscInt i, j; 74c4762a1bSJed Brown 75c4762a1bSJed Brown PetscFunctionBeginUser; 76c4762a1bSJed Brown /* Set values for independent variables and parameters */ 779566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 78c4762a1bSJed Brown x_a[0].setValue(x[0]); 79c4762a1bSJed Brown x_a[1].setValue(x[1]); 80c4762a1bSJed Brown mu_a.setValue(user->mu); 819566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 82c4762a1bSJed Brown 83c4762a1bSJed Brown /* Set seed matrix as 3x3 identity matrix */ 849371c9d4SSatish Balay x_a[0].setADValue(0, 1.); 859371c9d4SSatish Balay x_a[0].setADValue(1, 0.); 869371c9d4SSatish Balay x_a[0].setADValue(2, 0.); 879371c9d4SSatish Balay x_a[1].setADValue(0, 0.); 889371c9d4SSatish Balay x_a[1].setADValue(1, 1.); 899371c9d4SSatish Balay x_a[1].setADValue(2, 0.); 909371c9d4SSatish Balay mu_a.setADValue(0, 0.); 919371c9d4SSatish Balay mu_a.setADValue(1, 0.); 929371c9d4SSatish Balay mu_a.setADValue(2, 1.); 93c4762a1bSJed Brown 94c4762a1bSJed Brown /* Evaluate residual (on active variables) */ 959566063dSJacob Faibussowitsch PetscCall(EvaluateResidual(x_a, mu_a, f_a)); 96c4762a1bSJed Brown 97c4762a1bSJed Brown /* Extract derivatives */ 989566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(user->adctx->n, &J)); 99c4762a1bSJed Brown J[0] = (PetscScalar *)f_a[0].getADValue(); 100c4762a1bSJed Brown J[1] = (PetscScalar *)f_a[1].getADValue(); 101c4762a1bSJed Brown 102c4762a1bSJed Brown /* Set matrix values */ 103c4762a1bSJed Brown for (i = 0; i < user->adctx->m; i++) { 104*48a46eb9SPierre Jolivet for (j = 0; j < user->adctx->n; j++) PetscCall(MatSetValues(A, 1, &i, 1, &j, &J[i][j], INSERT_VALUES)); 105c4762a1bSJed Brown } 1069566063dSJacob Faibussowitsch PetscCall(PetscFree(J)); 1079566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 1089566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 109c4762a1bSJed Brown if (A != B) { 1109566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 1119566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 112c4762a1bSJed Brown } 113c4762a1bSJed Brown PetscFunctionReturn(0); 114c4762a1bSJed Brown } 115c4762a1bSJed Brown 116c4762a1bSJed Brown /* 117c4762a1bSJed Brown Compute the Jacobian w.r.t. mu using tapeless mode of ADOL-C. 118c4762a1bSJed Brown */ 1199371c9d4SSatish Balay static PetscErrorCode RHSJacobianP(TS ts, PetscReal t, Vec X, Mat A, void *ctx) { 120c4762a1bSJed Brown User user = (User)ctx; 121c4762a1bSJed Brown PetscScalar **J; 122c4762a1bSJed Brown PetscScalar *x; 123c4762a1bSJed Brown adouble f_a[2]; /* 'active' double for dependent variables */ 124c4762a1bSJed Brown adouble x_a[2], mu_a; /* 'active' doubles for independent variables */ 125c4762a1bSJed Brown PetscInt i, j = 0; 126c4762a1bSJed Brown 127c4762a1bSJed Brown PetscFunctionBeginUser; 128c4762a1bSJed Brown 129c4762a1bSJed Brown /* Set values for independent variables and parameters */ 1309566063dSJacob Faibussowitsch PetscCall(VecGetArray(X, &x)); 131c4762a1bSJed Brown x_a[0].setValue(x[0]); 132c4762a1bSJed Brown x_a[1].setValue(x[1]); 133c4762a1bSJed Brown mu_a.setValue(user->mu); 1349566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(X, &x)); 135c4762a1bSJed Brown 136c4762a1bSJed Brown /* Set seed matrix as 3x3 identity matrix */ 1379371c9d4SSatish Balay x_a[0].setADValue(0, 1.); 1389371c9d4SSatish Balay x_a[0].setADValue(1, 0.); 1399371c9d4SSatish Balay x_a[0].setADValue(2, 0.); 1409371c9d4SSatish Balay x_a[1].setADValue(0, 0.); 1419371c9d4SSatish Balay x_a[1].setADValue(1, 1.); 1429371c9d4SSatish Balay x_a[1].setADValue(2, 0.); 1439371c9d4SSatish Balay mu_a.setADValue(0, 0.); 1449371c9d4SSatish Balay mu_a.setADValue(1, 0.); 1459371c9d4SSatish Balay mu_a.setADValue(2, 1.); 146c4762a1bSJed Brown 147c4762a1bSJed Brown /* Evaluate residual (on active variables) */ 1489566063dSJacob Faibussowitsch PetscCall(EvaluateResidual(x_a, mu_a, f_a)); 149c4762a1bSJed Brown 150c4762a1bSJed Brown /* Extract derivatives */ 1519566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2, &J)); 152c4762a1bSJed Brown J[0] = (PetscScalar *)f_a[0].getADValue(); 153c4762a1bSJed Brown J[1] = (PetscScalar *)f_a[1].getADValue(); 154c4762a1bSJed Brown 155c4762a1bSJed Brown /* Set matrix values */ 156*48a46eb9SPierre Jolivet for (i = 0; i < user->adctx->m; i++) PetscCall(MatSetValues(A, 1, &i, 1, &j, &J[i][user->adctx->n], INSERT_VALUES)); 1579566063dSJacob Faibussowitsch PetscCall(PetscFree(J)); 1589566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 1599566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 160c4762a1bSJed Brown PetscFunctionReturn(0); 161c4762a1bSJed Brown } 162c4762a1bSJed Brown 163c4762a1bSJed Brown /* 164c4762a1bSJed Brown Monitor timesteps and use interpolation to output at integer multiples of 0.1 165c4762a1bSJed Brown */ 1669371c9d4SSatish Balay static PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal t, Vec X, void *ctx) { 167c4762a1bSJed Brown const PetscScalar *x; 168c4762a1bSJed Brown PetscReal tfinal, dt, tprev; 169c4762a1bSJed Brown User user = (User)ctx; 170c4762a1bSJed Brown 171c4762a1bSJed Brown PetscFunctionBeginUser; 1729566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts, &dt)); 1739566063dSJacob Faibussowitsch PetscCall(TSGetMaxTime(ts, &tfinal)); 1749566063dSJacob Faibussowitsch PetscCall(TSGetPrevTime(ts, &tprev)); 1759566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 17663a3b9bcSJacob 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]))); 1779566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "t %.6f (tprev = %.6f) \n", (double)t, (double)tprev)); 1789566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 179c4762a1bSJed Brown PetscFunctionReturn(0); 180c4762a1bSJed Brown } 181c4762a1bSJed Brown 1829371c9d4SSatish Balay int main(int argc, char **argv) { 183c4762a1bSJed Brown TS ts; /* nonlinear solver */ 184c4762a1bSJed Brown Vec x; /* solution, residual vectors */ 185c4762a1bSJed Brown Mat A; /* Jacobian matrix */ 186c4762a1bSJed Brown Mat Jacp; /* JacobianP matrix */ 187c4762a1bSJed Brown PetscInt steps; 188c4762a1bSJed Brown PetscReal ftime = 0.5; 189c4762a1bSJed Brown PetscBool monitor = PETSC_FALSE; 190c4762a1bSJed Brown PetscScalar *x_ptr; 191c4762a1bSJed Brown PetscMPIInt size; 192c4762a1bSJed Brown struct _n_User user; 193c4762a1bSJed Brown AdolcCtx *adctx; 194c4762a1bSJed Brown Vec lambda[2], mu[2]; 195c4762a1bSJed Brown 196c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 197c4762a1bSJed Brown Initialize program 198c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 199327415f7SBarry Smith PetscFunctionBeginUser; 2009566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 2019566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 20208401ef6SPierre Jolivet PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!"); 203c4762a1bSJed Brown 204c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 205c4762a1bSJed Brown Set runtime options and create AdolcCtx 206c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2079566063dSJacob Faibussowitsch PetscCall(PetscNew(&adctx)); 208c4762a1bSJed Brown user.mu = 1; 209c4762a1bSJed Brown user.next_output = 0.0; 2109371c9d4SSatish Balay adctx->m = 2; 2119371c9d4SSatish Balay adctx->n = 2; 2129371c9d4SSatish Balay adctx->p = 2; 213c4762a1bSJed Brown user.adctx = adctx; 214c4762a1bSJed Brown adtl::setNumDir(adctx->n + 1); /* #indep. variables, plus parameters */ 215c4762a1bSJed Brown 2169566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-mu", &user.mu, NULL)); 2179566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-monitor", &monitor, NULL)); 218c4762a1bSJed Brown 219c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 220c4762a1bSJed Brown Create necessary matrix and vectors, solve same ODE on every process 221c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2229566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 2239566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, 2, 2)); 2249566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 2259566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 2269566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &x, NULL)); 227c4762a1bSJed Brown 2289566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &Jacp)); 2299566063dSJacob Faibussowitsch PetscCall(MatSetSizes(Jacp, PETSC_DECIDE, PETSC_DECIDE, 2, 1)); 2309566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(Jacp)); 2319566063dSJacob Faibussowitsch PetscCall(MatSetUp(Jacp)); 232c4762a1bSJed Brown 233c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 234c4762a1bSJed Brown Create timestepping solver context 235c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2369566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 2379566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSRK)); 2389566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, NULL, RHSFunctionPassive, &user)); 239c4762a1bSJed Brown 240c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 241c4762a1bSJed Brown Set initial conditions 242c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2439566063dSJacob Faibussowitsch PetscCall(VecGetArray(x, &x_ptr)); 2449371c9d4SSatish Balay x_ptr[0] = 2; 2459371c9d4SSatish Balay x_ptr[1] = 0.66666654321; 2469566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x, &x_ptr)); 247c4762a1bSJed Brown 248c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 249c4762a1bSJed Brown Set RHS Jacobian for the adjoint integration 250c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2519566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(ts, A, A, RHSJacobian, &user)); 2529566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, ftime)); 2539566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP)); 254*48a46eb9SPierre Jolivet if (monitor) PetscCall(TSMonitorSet(ts, Monitor, &user, NULL)); 2559566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, .001)); 256c4762a1bSJed Brown 257c4762a1bSJed Brown /* 258c4762a1bSJed Brown Have the TS save its trajectory so that TSAdjointSolve() may be used 259c4762a1bSJed Brown */ 2609566063dSJacob Faibussowitsch PetscCall(TSSetSaveTrajectory(ts)); 261c4762a1bSJed Brown 262c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 263c4762a1bSJed Brown Set runtime options 264c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2659566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 266c4762a1bSJed Brown 267c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 268c4762a1bSJed Brown Solve nonlinear system 269c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2709566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, x)); 2719566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts, &ftime)); 2729566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &steps)); 27363a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "mu %g, steps %" PetscInt_FMT ", ftime %g\n", (double)user.mu, steps, (double)ftime)); 2749566063dSJacob Faibussowitsch PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD)); 275c4762a1bSJed Brown 276c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 277c4762a1bSJed Brown Start the Adjoint model 278c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2799566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &lambda[0], NULL)); 2809566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &lambda[1], NULL)); 281c4762a1bSJed Brown /* Reset initial conditions for the adjoint integration */ 2829566063dSJacob Faibussowitsch PetscCall(VecGetArray(lambda[0], &x_ptr)); 2839371c9d4SSatish Balay x_ptr[0] = 1.0; 2849371c9d4SSatish Balay x_ptr[1] = 0.0; 2859566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(lambda[0], &x_ptr)); 2869566063dSJacob Faibussowitsch PetscCall(VecGetArray(lambda[1], &x_ptr)); 2879371c9d4SSatish Balay x_ptr[0] = 0.0; 2889371c9d4SSatish Balay x_ptr[1] = 1.0; 2899566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(lambda[1], &x_ptr)); 290c4762a1bSJed Brown 2919566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(Jacp, &mu[0], NULL)); 2929566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(Jacp, &mu[1], NULL)); 2939566063dSJacob Faibussowitsch PetscCall(VecGetArray(mu[0], &x_ptr)); 294c4762a1bSJed Brown x_ptr[0] = 0.0; 2959566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(mu[0], &x_ptr)); 2969566063dSJacob Faibussowitsch PetscCall(VecGetArray(mu[1], &x_ptr)); 297c4762a1bSJed Brown x_ptr[0] = 0.0; 2989566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(mu[1], &x_ptr)); 2999566063dSJacob Faibussowitsch PetscCall(TSSetCostGradients(ts, 2, lambda, mu)); 300c4762a1bSJed Brown 301c4762a1bSJed Brown /* Set RHS JacobianP */ 3029566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobianP(ts, Jacp, RHSJacobianP, &user)); 303c4762a1bSJed Brown 3049566063dSJacob Faibussowitsch PetscCall(TSAdjointSolve(ts)); 305c4762a1bSJed Brown 3069566063dSJacob Faibussowitsch PetscCall(VecView(lambda[0], PETSC_VIEWER_STDOUT_WORLD)); 3079566063dSJacob Faibussowitsch PetscCall(VecView(lambda[1], PETSC_VIEWER_STDOUT_WORLD)); 3089566063dSJacob Faibussowitsch PetscCall(VecView(mu[0], PETSC_VIEWER_STDOUT_WORLD)); 3099566063dSJacob Faibussowitsch PetscCall(VecView(mu[1], PETSC_VIEWER_STDOUT_WORLD)); 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 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 3159566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 3169566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Jacp)); 3179566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 3189566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lambda[0])); 3199566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lambda[1])); 3209566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mu[0])); 3219566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mu[1])); 3229566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 3239566063dSJacob Faibussowitsch PetscCall(PetscFree(adctx)); 3249566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 325b122ec5aSJacob Faibussowitsch return 0; 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_tl_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_tl_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_tl_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_tl_4.out 352c4762a1bSJed Brown 353c4762a1bSJed Brown TEST*/ 354