1 2 static char help[] = "Solves the van der Pol equation and demonstrate IMEX.\n\ 3 Input parameters include:\n\ 4 -mu : stiffness parameter\n\n"; 5 6 /* ------------------------------------------------------------------------ 7 8 This program solves the van der Pol equation 9 y'' - \mu ((1-y^2)*y' - y) = 0 (1) 10 on the domain 0 <= x <= 1, with the boundary conditions 11 y(0) = 2, y'(0) = - 2/3 +10/(81*\mu) - 292/(2187*\mu^2), 12 This is a nonlinear equation. The well prepared initial condition gives errors that are not dominated by the first few steps of the method when \mu is large. 13 14 Notes: 15 This code demonstrates the TS solver interface to two variants of 16 linear problems, u_t = f(u,t), namely turning (1) into a system of 17 first order differential equations, 18 19 [ y' ] = [ z ] 20 [ z' ] [ \mu ((1 - y^2) z - y) ] 21 22 which then we can write as a vector equation 23 24 [ u_1' ] = [ u_2 ] (2) 25 [ u_2' ] [ \mu (1 - u_1^2) u_2 - u_1 ] 26 27 which is now in the desired form of u_t = f(u,t). One way that we 28 can split f(u,t) in (2) is to split by component, 29 30 [ u_1' ] = [ u_2 ] + [ 0 ] 31 [ u_2' ] [ 0 ] [ \mu ((1 - u_1^2) u_2 - u_1) ] 32 33 where 34 35 [ G(u,t) ] = [ u_2 ] 36 [ 0 ] 37 38 and 39 40 [ F(u',u,t) ] = [ u_1' ] - [ 0 ] 41 [ u_2' ] [ \mu ((1 - u_1^2) u_2 - u_1) ] 42 43 Using the definition of the Jacobian of F (from the PETSc user manual), 44 in the equation F(u',u,t) = G(u,t), 45 46 dF dF 47 J(F) = a * -- - -- 48 du' du 49 50 where d is the partial derivative. In this example, 51 52 dF [ 1 ; 0 ] 53 -- = [ ] 54 du' [ 0 ; 1 ] 55 56 dF [ 0 ; 0 ] 57 -- = [ ] 58 du [ -\mu (2*u_1*u_2 + 1); \mu (1 - u_1^2) ] 59 60 Hence, 61 62 [ a ; 0 ] 63 J(F) = [ ] 64 [ \mu (2*u_1*u_2 + 1); a - \mu (1 - u_1^2) ] 65 66 ------------------------------------------------------------------------- */ 67 68 #include <petscts.h> 69 70 typedef struct _n_User *User; 71 struct _n_User { 72 PetscReal mu; 73 PetscBool imex; 74 PetscReal next_output; 75 }; 76 77 /* 78 User-defined routines 79 */ 80 static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec X, Vec F, void *ctx) { 81 User user = (User)ctx; 82 PetscScalar *f; 83 const PetscScalar *x; 84 85 PetscFunctionBeginUser; 86 PetscCall(VecGetArrayRead(X, &x)); 87 PetscCall(VecGetArray(F, &f)); 88 f[0] = (user->imex ? x[1] : 0); 89 f[1] = 0.0; 90 PetscCall(VecRestoreArrayRead(X, &x)); 91 PetscCall(VecRestoreArray(F, &f)); 92 PetscFunctionReturn(0); 93 } 94 95 static PetscErrorCode IFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, void *ctx) { 96 User user = (User)ctx; 97 const PetscScalar *x, *xdot; 98 PetscScalar *f; 99 100 PetscFunctionBeginUser; 101 PetscCall(VecGetArrayRead(X, &x)); 102 PetscCall(VecGetArrayRead(Xdot, &xdot)); 103 PetscCall(VecGetArray(F, &f)); 104 f[0] = xdot[0] + (user->imex ? 0 : x[1]); 105 f[1] = xdot[1] - user->mu * ((1. - x[0] * x[0]) * x[1] - x[0]); 106 PetscCall(VecRestoreArrayRead(X, &x)); 107 PetscCall(VecRestoreArrayRead(Xdot, &xdot)); 108 PetscCall(VecRestoreArray(F, &f)); 109 PetscFunctionReturn(0); 110 } 111 112 static PetscErrorCode IJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat A, Mat B, void *ctx) { 113 User user = (User)ctx; 114 PetscReal mu = user->mu; 115 PetscInt rowcol[] = {0, 1}; 116 const PetscScalar *x; 117 PetscScalar J[2][2]; 118 119 PetscFunctionBeginUser; 120 PetscCall(VecGetArrayRead(X, &x)); 121 J[0][0] = a; 122 J[0][1] = (user->imex ? 0 : 1.); 123 J[1][0] = mu * (2. * x[0] * x[1] + 1.); 124 J[1][1] = a - mu * (1. - x[0] * x[0]); 125 PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES)); 126 PetscCall(VecRestoreArrayRead(X, &x)); 127 128 PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 129 PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 130 if (A != B) { 131 PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 132 PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 133 } 134 PetscFunctionReturn(0); 135 } 136 137 static PetscErrorCode RegisterMyARK2(void) { 138 PetscFunctionBeginUser; 139 { 140 const PetscReal A[3][3] = 141 { 142 {0, 0, 0}, 143 {0.41421356237309504880, 0, 0}, 144 {0.75, 0.25, 0} 145 }, 146 At[3][3] = {{0, 0, 0}, {0.12132034355964257320, 0.29289321881345247560, 0}, {0.20710678118654752440, 0.50000000000000000000, 0.29289321881345247560}}, *bembedt = NULL, *bembed = NULL; 147 PetscCall(TSARKIMEXRegister("myark2", 2, 3, &At[0][0], NULL, NULL, &A[0][0], NULL, NULL, bembedt, bembed, 0, NULL, NULL)); 148 } 149 PetscFunctionReturn(0); 150 } 151 152 /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */ 153 static PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal t, Vec X, void *ctx) { 154 const PetscScalar *x; 155 PetscReal tfinal, dt; 156 User user = (User)ctx; 157 Vec interpolatedX; 158 159 PetscFunctionBeginUser; 160 PetscCall(TSGetTimeStep(ts, &dt)); 161 PetscCall(TSGetMaxTime(ts, &tfinal)); 162 163 while (user->next_output <= t && user->next_output <= tfinal) { 164 PetscCall(VecDuplicate(X, &interpolatedX)); 165 PetscCall(TSInterpolate(ts, user->next_output, interpolatedX)); 166 PetscCall(VecGetArrayRead(interpolatedX, &x)); 167 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]))); 168 PetscCall(VecRestoreArrayRead(interpolatedX, &x)); 169 PetscCall(VecDestroy(&interpolatedX)); 170 171 user->next_output += 0.1; 172 } 173 PetscFunctionReturn(0); 174 } 175 176 int main(int argc, char **argv) { 177 TS ts; /* nonlinear solver */ 178 Vec x; /* solution, residual vectors */ 179 Mat A; /* Jacobian matrix */ 180 PetscInt steps; 181 PetscReal ftime = 0.5; 182 PetscBool monitor = PETSC_FALSE; 183 PetscScalar *x_ptr; 184 PetscMPIInt size; 185 struct _n_User user; 186 187 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 188 Initialize program 189 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 190 PetscFunctionBeginUser; 191 PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 192 PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 193 PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!"); 194 195 PetscCall(RegisterMyARK2()); 196 197 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 198 Set runtime options 199 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 200 user.mu = 1000.0; 201 user.imex = PETSC_TRUE; 202 user.next_output = 0.0; 203 204 PetscCall(PetscOptionsGetReal(NULL, NULL, "-mu", &user.mu, NULL)); 205 PetscCall(PetscOptionsGetBool(NULL, NULL, "-imex", &user.imex, NULL)); 206 PetscCall(PetscOptionsGetBool(NULL, NULL, "-monitor", &monitor, NULL)); 207 208 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 209 Create necessary matrix and vectors, solve same ODE on every process 210 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 211 PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 212 PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, 2, 2)); 213 PetscCall(MatSetFromOptions(A)); 214 PetscCall(MatSetUp(A)); 215 PetscCall(MatCreateVecs(A, &x, NULL)); 216 217 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 218 Create timestepping solver context 219 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 220 PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 221 PetscCall(TSSetType(ts, TSBEULER)); 222 PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, &user)); 223 PetscCall(TSSetIFunction(ts, NULL, IFunction, &user)); 224 PetscCall(TSSetIJacobian(ts, A, A, IJacobian, &user)); 225 PetscCall(TSSetMaxTime(ts, ftime)); 226 PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER)); 227 if (monitor) PetscCall(TSMonitorSet(ts, Monitor, &user, NULL)); 228 229 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 230 Set initial conditions 231 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 232 PetscCall(VecGetArray(x, &x_ptr)); 233 x_ptr[0] = 2.0; 234 x_ptr[1] = -2.0 / 3.0 + 10.0 / (81.0 * user.mu) - 292.0 / (2187.0 * user.mu * user.mu); 235 PetscCall(VecRestoreArray(x, &x_ptr)); 236 PetscCall(TSSetTimeStep(ts, 0.01)); 237 238 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 239 Set runtime options 240 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 241 PetscCall(TSSetFromOptions(ts)); 242 243 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 244 Solve nonlinear system 245 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 246 PetscCall(TSSolve(ts, x)); 247 PetscCall(TSGetSolveTime(ts, &ftime)); 248 PetscCall(TSGetStepNumber(ts, &steps)); 249 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "mu %g, steps %" PetscInt_FMT ", ftime %g\n", (double)user.mu, steps, (double)ftime)); 250 PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD)); 251 252 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 253 Free work space. All PETSc objects should be destroyed when they 254 are no longer needed. 255 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 256 PetscCall(MatDestroy(&A)); 257 PetscCall(VecDestroy(&x)); 258 PetscCall(TSDestroy(&ts)); 259 260 PetscCall(PetscFinalize()); 261 return 0; 262 } 263 264 /*TEST 265 266 test: 267 args: -ts_type arkimex -ts_arkimex_type myark2 -ts_adapt_type none 268 requires: !single 269 270 TEST*/ 271