1 2 static char help[] = "Solves the van der Pol equation.\n\ 3 Input parameters include:\n"; 4 5 /* ------------------------------------------------------------------------ 6 7 This program solves the van der Pol DAE ODE equivalent 8 y' = z (1) 9 z' = \mu ((1-y^2)z-y) 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 and 13 \mu = 10^6 ( y'(0) ~ -0.6666665432100101). 14 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. 15 16 Notes: 17 This code demonstrates the TS solver interface to an ODE -- RHSFunction for explicit form and IFunction for implicit form. 18 19 ------------------------------------------------------------------------- */ 20 21 #include <petscts.h> 22 23 typedef struct _n_User *User; 24 struct _n_User { 25 PetscReal mu; 26 PetscReal next_output; 27 }; 28 29 /* 30 User-defined routines 31 */ 32 static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec X, Vec F, void *ctx) { 33 User user = (User)ctx; 34 PetscScalar *f; 35 const PetscScalar *x; 36 37 PetscFunctionBeginUser; 38 PetscCall(VecGetArrayRead(X, &x)); 39 PetscCall(VecGetArray(F, &f)); 40 f[0] = x[1]; 41 f[1] = user->mu * (1. - x[0] * x[0]) * x[1] - x[0]; 42 PetscCall(VecRestoreArrayRead(X, &x)); 43 PetscCall(VecRestoreArray(F, &f)); 44 PetscFunctionReturn(0); 45 } 46 47 static PetscErrorCode IFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, void *ctx) { 48 User user = (User)ctx; 49 const PetscScalar *x, *xdot; 50 PetscScalar *f; 51 52 PetscFunctionBeginUser; 53 PetscCall(VecGetArrayRead(X, &x)); 54 PetscCall(VecGetArrayRead(Xdot, &xdot)); 55 PetscCall(VecGetArray(F, &f)); 56 f[0] = xdot[0] - x[1]; 57 f[1] = xdot[1] - user->mu * ((1.0 - x[0] * x[0]) * x[1] - x[0]); 58 PetscCall(VecRestoreArrayRead(X, &x)); 59 PetscCall(VecRestoreArrayRead(Xdot, &xdot)); 60 PetscCall(VecRestoreArray(F, &f)); 61 PetscFunctionReturn(0); 62 } 63 64 static PetscErrorCode IJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat A, Mat B, void *ctx) { 65 User user = (User)ctx; 66 PetscInt rowcol[] = {0, 1}; 67 const PetscScalar *x; 68 PetscScalar J[2][2]; 69 70 PetscFunctionBeginUser; 71 PetscCall(VecGetArrayRead(X, &x)); 72 J[0][0] = a; 73 J[0][1] = -1.0; 74 J[1][0] = user->mu * (2.0 * x[0] * x[1] + 1.0); 75 J[1][1] = a - user->mu * (1.0 - x[0] * x[0]); 76 PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES)); 77 PetscCall(VecRestoreArrayRead(X, &x)); 78 79 PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 80 PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 81 if (A != B) { 82 PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 83 PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 84 } 85 PetscFunctionReturn(0); 86 } 87 88 /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */ 89 static PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal t, Vec X, void *ctx) { 90 const PetscScalar *x; 91 PetscReal tfinal, dt; 92 User user = (User)ctx; 93 Vec interpolatedX; 94 95 PetscFunctionBeginUser; 96 PetscCall(TSGetTimeStep(ts, &dt)); 97 PetscCall(TSGetMaxTime(ts, &tfinal)); 98 99 while (user->next_output <= t && user->next_output <= tfinal) { 100 PetscCall(VecDuplicate(X, &interpolatedX)); 101 PetscCall(TSInterpolate(ts, user->next_output, interpolatedX)); 102 PetscCall(VecGetArrayRead(interpolatedX, &x)); 103 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]))); 104 PetscCall(VecRestoreArrayRead(interpolatedX, &x)); 105 PetscCall(VecDestroy(&interpolatedX)); 106 user->next_output += 0.1; 107 } 108 PetscFunctionReturn(0); 109 } 110 111 int main(int argc, char **argv) { 112 TS ts; /* nonlinear solver */ 113 Vec x; /* solution, residual vectors */ 114 Mat A; /* Jacobian matrix */ 115 PetscInt steps; 116 PetscReal ftime = 0.5; 117 PetscBool monitor = PETSC_FALSE, implicitform = PETSC_TRUE; 118 PetscScalar *x_ptr; 119 PetscMPIInt size; 120 struct _n_User user; 121 122 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 123 Initialize program 124 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 125 PetscFunctionBeginUser; 126 PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 127 PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 128 PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!"); 129 130 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 131 Set runtime options 132 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 133 user.next_output = 0.0; 134 user.mu = 1.0e3; 135 PetscCall(PetscOptionsGetBool(NULL, NULL, "-monitor", &monitor, NULL)); 136 PetscCall(PetscOptionsGetBool(NULL, NULL, "-implicitform", &implicitform, NULL)); 137 PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Physical parameters", NULL); 138 PetscCall(PetscOptionsReal("-mu", "Stiffness parameter", "<1.0e6>", user.mu, &user.mu, NULL)); 139 PetscOptionsEnd(); 140 141 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 142 Create necessary matrix and vectors, solve same ODE on every process 143 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 144 PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 145 PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, 2, 2)); 146 PetscCall(MatSetFromOptions(A)); 147 PetscCall(MatSetUp(A)); 148 149 PetscCall(MatCreateVecs(A, &x, NULL)); 150 151 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 152 Create timestepping solver context 153 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 154 PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 155 if (implicitform) { 156 PetscCall(TSSetIFunction(ts, NULL, IFunction, &user)); 157 PetscCall(TSSetIJacobian(ts, A, A, IJacobian, &user)); 158 PetscCall(TSSetType(ts, TSBEULER)); 159 } else { 160 PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, &user)); 161 PetscCall(TSSetType(ts, TSRK)); 162 } 163 PetscCall(TSSetMaxTime(ts, ftime)); 164 PetscCall(TSSetTimeStep(ts, 0.001)); 165 PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER)); 166 if (monitor) PetscCall(TSMonitorSet(ts, Monitor, &user, NULL)); 167 168 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 169 Set initial conditions 170 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 171 PetscCall(VecGetArray(x, &x_ptr)); 172 x_ptr[0] = 2.0; 173 x_ptr[1] = -2.0 / 3.0 + 10.0 / (81.0 * user.mu) - 292.0 / (2187.0 * user.mu * user.mu); 174 PetscCall(VecRestoreArray(x, &x_ptr)); 175 176 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 177 Set runtime options 178 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 179 PetscCall(TSSetFromOptions(ts)); 180 181 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 182 Solve nonlinear system 183 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 184 PetscCall(TSSolve(ts, x)); 185 PetscCall(TSGetSolveTime(ts, &ftime)); 186 PetscCall(TSGetStepNumber(ts, &steps)); 187 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "steps %" PetscInt_FMT ", ftime %g\n", steps, (double)ftime)); 188 PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD)); 189 190 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 191 Free work space. All PETSc objects should be destroyed when they 192 are no longer needed. 193 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 194 PetscCall(MatDestroy(&A)); 195 PetscCall(VecDestroy(&x)); 196 PetscCall(TSDestroy(&ts)); 197 198 PetscCall(PetscFinalize()); 199 return (0); 200 } 201 202 /*TEST 203 204 test: 205 requires: !single 206 args: -mu 1e6 207 208 test: 209 requires: !single 210 suffix: 2 211 args: -implicitform false -ts_type rk -ts_rk_type 5dp -ts_adapt_type dsp 212 213 test: 214 requires: !single 215 suffix: 3 216 args: -implicitform false -ts_type rk -ts_rk_type 5dp -ts_adapt_type dsp -ts_adapt_dsp_filter H0312 217 218 TEST*/ 219