1 static char help[] = "Performs adjoint sensitivity analysis for the van der Pol equation.\n\ 2 Input parameters include:\n\ 3 -mu : stiffness parameter\n\n"; 4 5 /* ------------------------------------------------------------------------ 6 7 This program solves the van der Pol equation 8 y'' - \mu (1-y^2)*y' + y = 0 (1) 9 on the domain 0 <= x <= 1, with the boundary conditions 10 y(0) = 2, y'(0) = 0, 11 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. 12 13 Notes: 14 This code demonstrates the TSForward interface to a system of ordinary differential equations (ODEs) in the form of u_t = f(u,t). 15 16 (1) can be turned into a system of first order ODEs 17 [ y' ] = [ z ] 18 [ z' ] [ \mu (1 - y^2) z - y ] 19 20 which then we can write as a vector equation 21 22 [ u_1' ] = [ u_2 ] (2) 23 [ u_2' ] [ \mu (1 - u_1^2) u_2 - u_1 ] 24 25 which is now in the form of u_t = F(u,t). 26 27 The user provides the right-hand-side function 28 29 [ f(u,t) ] = [ u_2 ] 30 [ \mu (1 - u_1^2) u_2 - u_1 ] 31 32 the Jacobian function 33 34 df [ 0 ; 1 ] 35 -- = [ ] 36 du [ -2 \mu u_1*u_2 - 1; \mu (1 - u_1^2) ] 37 38 and the JacobainP (the Jacobian w.r.t. parameter) function 39 40 df [ 0; 0; 0 ] 41 --- = [ ] 42 d\mu [ 0; 0; (1 - u_1^2) u_2 ] 43 44 ------------------------------------------------------------------------- */ 45 46 #include <petscts.h> 47 #include <petscmat.h> 48 typedef struct _n_User *User; 49 struct _n_User { 50 PetscReal mu; 51 PetscReal next_output; 52 PetscReal tprev; 53 }; 54 55 /* 56 User-defined routines 57 */ 58 static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec X, Vec F, void *ctx) { 59 User user = (User)ctx; 60 PetscScalar *f; 61 const PetscScalar *x; 62 63 PetscFunctionBeginUser; 64 PetscCall(VecGetArrayRead(X, &x)); 65 PetscCall(VecGetArray(F, &f)); 66 f[0] = x[1]; 67 f[1] = user->mu * (1. - x[0] * x[0]) * x[1] - x[0]; 68 PetscCall(VecRestoreArrayRead(X, &x)); 69 PetscCall(VecRestoreArray(F, &f)); 70 PetscFunctionReturn(0); 71 } 72 73 static PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec X, Mat A, Mat B, void *ctx) { 74 User user = (User)ctx; 75 PetscReal mu = user->mu; 76 PetscInt rowcol[] = {0, 1}; 77 PetscScalar J[2][2]; 78 const PetscScalar *x; 79 80 PetscFunctionBeginUser; 81 PetscCall(VecGetArrayRead(X, &x)); 82 J[0][0] = 0; 83 J[1][0] = -2. * mu * x[1] * x[0] - 1.; 84 J[0][1] = 1.0; 85 J[1][1] = mu * (1.0 - x[0] * x[0]); 86 PetscCall(MatSetValues(A, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES)); 87 PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 88 PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 89 if (A != B) { 90 PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 91 PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 92 } 93 PetscCall(VecRestoreArrayRead(X, &x)); 94 PetscFunctionReturn(0); 95 } 96 97 static PetscErrorCode RHSJacobianP(TS ts, PetscReal t, Vec X, Mat A, void *ctx) { 98 PetscInt row[] = {0, 1}, col[] = {2}; 99 PetscScalar J[2][1]; 100 const PetscScalar *x; 101 102 PetscFunctionBeginUser; 103 PetscCall(VecGetArrayRead(X, &x)); 104 J[0][0] = 0; 105 J[1][0] = (1. - x[0] * x[0]) * x[1]; 106 PetscCall(VecRestoreArrayRead(X, &x)); 107 PetscCall(MatSetValues(A, 2, row, 1, col, &J[0][0], INSERT_VALUES)); 108 109 PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 110 PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 111 PetscFunctionReturn(0); 112 } 113 114 /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */ 115 static PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal t, Vec X, void *ctx) { 116 const PetscScalar *x; 117 PetscReal tfinal, dt, tprev; 118 User user = (User)ctx; 119 120 PetscFunctionBeginUser; 121 PetscCall(TSGetTimeStep(ts, &dt)); 122 PetscCall(TSGetMaxTime(ts, &tfinal)); 123 PetscCall(TSGetPrevTime(ts, &tprev)); 124 PetscCall(VecGetArrayRead(X, &x)); 125 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]))); 126 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "t %.6f (tprev = %.6f) \n", (double)t, (double)tprev)); 127 PetscCall(VecRestoreArrayRead(X, &x)); 128 PetscFunctionReturn(0); 129 } 130 131 int main(int argc, char **argv) { 132 TS ts; /* nonlinear solver */ 133 Vec x; /* solution, residual vectors */ 134 Mat A; /* Jacobian matrix */ 135 Mat Jacp; /* JacobianP matrix */ 136 PetscInt steps; 137 PetscReal ftime = 0.5; 138 PetscBool monitor = PETSC_FALSE; 139 PetscScalar *x_ptr; 140 PetscMPIInt size; 141 struct _n_User user; 142 Mat sp; 143 144 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 145 Initialize program 146 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 147 PetscFunctionBeginUser; 148 PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 149 PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 150 PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!"); 151 152 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 153 Set runtime options 154 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 155 user.mu = 1; 156 user.next_output = 0.0; 157 158 PetscCall(PetscOptionsGetReal(NULL, NULL, "-mu", &user.mu, NULL)); 159 PetscCall(PetscOptionsGetBool(NULL, NULL, "-monitor", &monitor, NULL)); 160 161 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 162 Create necessary matrix and vectors, solve same ODE on every process 163 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 164 PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 165 PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, 2, 2)); 166 PetscCall(MatSetFromOptions(A)); 167 PetscCall(MatSetUp(A)); 168 PetscCall(MatCreateVecs(A, &x, NULL)); 169 170 PetscCall(MatCreate(PETSC_COMM_WORLD, &Jacp)); 171 PetscCall(MatSetSizes(Jacp, PETSC_DECIDE, PETSC_DECIDE, 2, 3)); 172 PetscCall(MatSetFromOptions(Jacp)); 173 PetscCall(MatSetUp(Jacp)); 174 175 PetscCall(MatCreateDense(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, 2, 3, NULL, &sp)); 176 PetscCall(MatZeroEntries(sp)); 177 PetscCall(MatShift(sp, 1.0)); 178 179 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 180 Create timestepping solver context 181 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 182 PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 183 PetscCall(TSSetType(ts, TSRK)); 184 PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, &user)); 185 /* Set RHS Jacobian for the adjoint integration */ 186 PetscCall(TSSetRHSJacobian(ts, A, A, RHSJacobian, &user)); 187 PetscCall(TSSetMaxTime(ts, ftime)); 188 PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP)); 189 if (monitor) { PetscCall(TSMonitorSet(ts, Monitor, &user, NULL)); } 190 PetscCall(TSForwardSetSensitivities(ts, 3, sp)); 191 PetscCall(TSSetRHSJacobianP(ts, Jacp, RHSJacobianP, &user)); 192 193 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 194 Set initial conditions 195 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 196 PetscCall(VecGetArray(x, &x_ptr)); 197 198 x_ptr[0] = 2; 199 x_ptr[1] = 0.66666654321; 200 PetscCall(VecRestoreArray(x, &x_ptr)); 201 PetscCall(TSSetTimeStep(ts, .001)); 202 203 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 204 Set runtime options 205 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 206 PetscCall(TSSetFromOptions(ts)); 207 208 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 209 Solve nonlinear system 210 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 211 PetscCall(TSSolve(ts, x)); 212 PetscCall(TSGetSolveTime(ts, &ftime)); 213 PetscCall(TSGetStepNumber(ts, &steps)); 214 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "mu %g, steps %" PetscInt_FMT ", ftime %g\n", (double)user.mu, steps, (double)ftime)); 215 PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD)); 216 217 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n forward sensitivity: d[y(tf) z(tf)]/d[y0 z0 mu]\n")); 218 PetscCall(MatView(sp, PETSC_VIEWER_STDOUT_WORLD)); 219 220 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 221 Free work space. All PETSc objects should be destroyed when they 222 are no longer needed. 223 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 224 PetscCall(MatDestroy(&A)); 225 PetscCall(MatDestroy(&Jacp)); 226 PetscCall(VecDestroy(&x)); 227 PetscCall(MatDestroy(&sp)); 228 PetscCall(TSDestroy(&ts)); 229 PetscCall(PetscFinalize()); 230 return 0; 231 } 232 233 /*TEST 234 235 test: 236 args: -monitor 0 -ts_adapt_type none 237 238 TEST*/ 239