1c4762a1bSJed Brown static char help[] = "Performs adjoint sensitivity analysis for 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 7c4762a1bSJed Brown This program solves the van der Pol equation 8c4762a1bSJed Brown y'' - \mu (1-y^2)*y' + y = 0 (1) 9c4762a1bSJed Brown on the domain 0 <= x <= 1, with the boundary conditions 10c4762a1bSJed Brown y(0) = 2, y'(0) = 0, 11c4762a1bSJed Brown and computes the sensitivities of the final solution w.r.t. initial conditions and parameter \mu with an explicit Runge-Kutta method and its discrete tangent linear model. 12c4762a1bSJed Brown 13c4762a1bSJed Brown Notes: 145ab1ac2bSHong Zhang This code demonstrates the TSForward interface to a system of ordinary differential equations (ODEs) in the form of u_t = f(u,t). 15c4762a1bSJed Brown 16c4762a1bSJed Brown (1) can be turned into a system of first order ODEs 17c4762a1bSJed Brown [ y' ] = [ z ] 18c4762a1bSJed Brown [ z' ] [ \mu (1 - y^2) z - y ] 19c4762a1bSJed Brown 20c4762a1bSJed Brown which then we can write as a vector equation 21c4762a1bSJed Brown 22c4762a1bSJed Brown [ u_1' ] = [ u_2 ] (2) 23c4762a1bSJed Brown [ u_2' ] [ \mu (1 - u_1^2) u_2 - u_1 ] 24c4762a1bSJed Brown 25c4762a1bSJed Brown which is now in the form of u_t = F(u,t). 26c4762a1bSJed Brown 27c4762a1bSJed Brown The user provides the right-hand-side function 28c4762a1bSJed Brown 295ab1ac2bSHong Zhang [ f(u,t) ] = [ u_2 ] 30c4762a1bSJed Brown [ \mu (1 - u_1^2) u_2 - u_1 ] 31c4762a1bSJed Brown 32c4762a1bSJed Brown the Jacobian function 33c4762a1bSJed Brown 345ab1ac2bSHong Zhang df [ 0 ; 1 ] 35c4762a1bSJed Brown -- = [ ] 36c4762a1bSJed Brown du [ -2 \mu u_1*u_2 - 1; \mu (1 - u_1^2) ] 37c4762a1bSJed Brown 38c4762a1bSJed Brown and the JacobainP (the Jacobian w.r.t. parameter) function 39c4762a1bSJed Brown 405ab1ac2bSHong Zhang df [ 0; 0; 0 ] 41c4762a1bSJed Brown --- = [ ] 42c4762a1bSJed Brown d\mu [ 0; 0; (1 - u_1^2) u_2 ] 43c4762a1bSJed Brown 44c4762a1bSJed Brown ------------------------------------------------------------------------- */ 45c4762a1bSJed Brown 46c4762a1bSJed Brown #include <petscts.h> 47c4762a1bSJed Brown #include <petscmat.h> 48c4762a1bSJed Brown typedef struct _n_User *User; 49c4762a1bSJed Brown struct _n_User { 50c4762a1bSJed Brown PetscReal mu; 51c4762a1bSJed Brown PetscReal next_output; 52c4762a1bSJed Brown PetscReal tprev; 53c4762a1bSJed Brown }; 54c4762a1bSJed Brown 55c4762a1bSJed Brown /* 560e3d61c9SBarry Smith User-defined routines 57c4762a1bSJed Brown */ 589371c9d4SSatish Balay static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec X, Vec F, void *ctx) { 59c4762a1bSJed Brown User user = (User)ctx; 60c4762a1bSJed Brown PetscScalar *f; 61c4762a1bSJed Brown const PetscScalar *x; 62c4762a1bSJed Brown 63c4762a1bSJed Brown PetscFunctionBeginUser; 649566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 659566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 66c4762a1bSJed Brown f[0] = x[1]; 67c4762a1bSJed Brown f[1] = user->mu * (1. - x[0] * x[0]) * x[1] - x[0]; 689566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 699566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 70c4762a1bSJed Brown PetscFunctionReturn(0); 71c4762a1bSJed Brown } 72c4762a1bSJed Brown 739371c9d4SSatish Balay static PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec X, Mat A, Mat B, void *ctx) { 74c4762a1bSJed Brown User user = (User)ctx; 75c4762a1bSJed Brown PetscReal mu = user->mu; 76c4762a1bSJed Brown PetscInt rowcol[] = {0, 1}; 77c4762a1bSJed Brown PetscScalar J[2][2]; 78c4762a1bSJed Brown const PetscScalar *x; 79c4762a1bSJed Brown 80c4762a1bSJed Brown PetscFunctionBeginUser; 819566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 82c4762a1bSJed Brown J[0][0] = 0; 83c4762a1bSJed Brown J[1][0] = -2. * mu * x[1] * x[0] - 1.; 84c4762a1bSJed Brown J[0][1] = 1.0; 85c4762a1bSJed Brown J[1][1] = mu * (1.0 - x[0] * x[0]); 869566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES)); 879566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 889566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 89c4762a1bSJed Brown if (A != B) { 909566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 919566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 92c4762a1bSJed Brown } 939566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 94c4762a1bSJed Brown PetscFunctionReturn(0); 95c4762a1bSJed Brown } 96c4762a1bSJed Brown 979371c9d4SSatish Balay static PetscErrorCode RHSJacobianP(TS ts, PetscReal t, Vec X, Mat A, void *ctx) { 98c4762a1bSJed Brown PetscInt row[] = {0, 1}, col[] = {2}; 99c4762a1bSJed Brown PetscScalar J[2][1]; 100c4762a1bSJed Brown const PetscScalar *x; 101c4762a1bSJed Brown 102c4762a1bSJed Brown PetscFunctionBeginUser; 1039566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 104c4762a1bSJed Brown J[0][0] = 0; 105c4762a1bSJed Brown J[1][0] = (1. - x[0] * x[0]) * x[1]; 1069566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 1079566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 2, row, 1, col, &J[0][0], INSERT_VALUES)); 108c4762a1bSJed Brown 1099566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 1109566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 111c4762a1bSJed Brown PetscFunctionReturn(0); 112c4762a1bSJed Brown } 113c4762a1bSJed Brown 114c4762a1bSJed Brown /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */ 1159371c9d4SSatish Balay static PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal t, Vec X, void *ctx) { 116c4762a1bSJed Brown const PetscScalar *x; 117c4762a1bSJed Brown PetscReal tfinal, dt, tprev; 118c4762a1bSJed Brown User user = (User)ctx; 119c4762a1bSJed Brown 120c4762a1bSJed Brown PetscFunctionBeginUser; 1219566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts, &dt)); 1229566063dSJacob Faibussowitsch PetscCall(TSGetMaxTime(ts, &tfinal)); 1239566063dSJacob Faibussowitsch PetscCall(TSGetPrevTime(ts, &tprev)); 1249566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 12563a3b9bcSJacob 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]))); 1269566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "t %.6f (tprev = %.6f) \n", (double)t, (double)tprev)); 1279566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 128c4762a1bSJed Brown PetscFunctionReturn(0); 129c4762a1bSJed Brown } 130c4762a1bSJed Brown 1319371c9d4SSatish Balay int main(int argc, char **argv) { 132c4762a1bSJed Brown TS ts; /* nonlinear solver */ 133c4762a1bSJed Brown Vec x; /* solution, residual vectors */ 134c4762a1bSJed Brown Mat A; /* Jacobian matrix */ 135c4762a1bSJed Brown Mat Jacp; /* JacobianP matrix */ 136c4762a1bSJed Brown PetscInt steps; 137c4762a1bSJed Brown PetscReal ftime = 0.5; 138c4762a1bSJed Brown PetscBool monitor = PETSC_FALSE; 139c4762a1bSJed Brown PetscScalar *x_ptr; 140c4762a1bSJed Brown PetscMPIInt size; 141c4762a1bSJed Brown struct _n_User user; 142c4762a1bSJed Brown Mat sp; 143c4762a1bSJed Brown 144c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 145c4762a1bSJed Brown Initialize program 146c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 147327415f7SBarry Smith PetscFunctionBeginUser; 1489566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 1499566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 1503c633725SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!"); 151c4762a1bSJed Brown 152c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 153c4762a1bSJed Brown Set runtime options 154c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 155c4762a1bSJed Brown user.mu = 1; 156c4762a1bSJed Brown user.next_output = 0.0; 157c4762a1bSJed Brown 1589566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-mu", &user.mu, NULL)); 1599566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-monitor", &monitor, NULL)); 160c4762a1bSJed Brown 161c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 162c4762a1bSJed Brown Create necessary matrix and vectors, solve same ODE on every process 163c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1649566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 1659566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, 2, 2)); 1669566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 1679566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 1689566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &x, NULL)); 169c4762a1bSJed Brown 1709566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &Jacp)); 1719566063dSJacob Faibussowitsch PetscCall(MatSetSizes(Jacp, PETSC_DECIDE, PETSC_DECIDE, 2, 3)); 1729566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(Jacp)); 1739566063dSJacob Faibussowitsch PetscCall(MatSetUp(Jacp)); 174c4762a1bSJed Brown 1759566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, 2, 3, NULL, &sp)); 1769566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(sp)); 1779566063dSJacob Faibussowitsch PetscCall(MatShift(sp, 1.0)); 178c4762a1bSJed Brown 179c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 180c4762a1bSJed Brown Create timestepping solver context 181c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1829566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 1839566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSRK)); 1849566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, &user)); 185c4762a1bSJed Brown /* Set RHS Jacobian for the adjoint integration */ 1869566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(ts, A, A, RHSJacobian, &user)); 1879566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, ftime)); 1889566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP)); 189*48a46eb9SPierre Jolivet if (monitor) PetscCall(TSMonitorSet(ts, Monitor, &user, NULL)); 1909566063dSJacob Faibussowitsch PetscCall(TSForwardSetSensitivities(ts, 3, sp)); 1919566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobianP(ts, Jacp, RHSJacobianP, &user)); 192c4762a1bSJed Brown 193c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 194c4762a1bSJed Brown Set initial conditions 195c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1969566063dSJacob Faibussowitsch PetscCall(VecGetArray(x, &x_ptr)); 197c4762a1bSJed Brown 1989371c9d4SSatish Balay x_ptr[0] = 2; 1999371c9d4SSatish Balay x_ptr[1] = 0.66666654321; 2009566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x, &x_ptr)); 2019566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, .001)); 202c4762a1bSJed Brown 203c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 204c4762a1bSJed Brown Set runtime options 205c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2069566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 207c4762a1bSJed Brown 208c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 209c4762a1bSJed Brown Solve nonlinear system 210c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2119566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, x)); 2129566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts, &ftime)); 2139566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &steps)); 21463a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "mu %g, steps %" PetscInt_FMT ", ftime %g\n", (double)user.mu, steps, (double)ftime)); 2159566063dSJacob Faibussowitsch PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD)); 216c4762a1bSJed Brown 2179566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n forward sensitivity: d[y(tf) z(tf)]/d[y0 z0 mu]\n")); 2189566063dSJacob Faibussowitsch PetscCall(MatView(sp, PETSC_VIEWER_STDOUT_WORLD)); 219c4762a1bSJed Brown 220c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 221c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 222c4762a1bSJed Brown are no longer needed. 223c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2249566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 2259566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Jacp)); 2269566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 2279566063dSJacob Faibussowitsch PetscCall(MatDestroy(&sp)); 2289566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 2299566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 230b122ec5aSJacob Faibussowitsch return 0; 231c4762a1bSJed Brown } 232c4762a1bSJed Brown 233c4762a1bSJed Brown /*TEST 234c4762a1bSJed Brown 235c4762a1bSJed Brown test: 236c4762a1bSJed Brown args: -monitor 0 -ts_adapt_type none 237c4762a1bSJed Brown 238c4762a1bSJed Brown TEST*/ 239