1c4762a1bSJed Brown static char help[] = "Solves the van der Pol equation.\n\ 2c4762a1bSJed Brown Input parameters include:\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown /* ------------------------------------------------------------------------ 5c4762a1bSJed Brown 6c4762a1bSJed Brown This program solves the van der Pol DAE ODE equivalent 7c4762a1bSJed Brown y' = z (1) 8c4762a1bSJed Brown z' = \mu ((1-y^2)z-y) 9c4762a1bSJed Brown on the domain 0 <= x <= 1, with the boundary conditions 10c4762a1bSJed Brown y(0) = 2, y'(0) = - 2/3 +10/(81*\mu) - 292/(2187*\mu^2), 11c4762a1bSJed Brown and 12c4762a1bSJed Brown \mu = 10^6 ( y'(0) ~ -0.6666665432100101). 13c4762a1bSJed Brown 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. 14c4762a1bSJed Brown 15c4762a1bSJed Brown Notes: 16c4762a1bSJed Brown This code demonstrates the TS solver interface to an ODE -- RHSFunction for explicit form and IFunction for implicit form. 17c4762a1bSJed Brown 18c4762a1bSJed Brown ------------------------------------------------------------------------- */ 19c4762a1bSJed Brown 20c4762a1bSJed Brown #include <petscts.h> 21c4762a1bSJed Brown 22c4762a1bSJed Brown typedef struct _n_User *User; 23c4762a1bSJed Brown struct _n_User { 24c4762a1bSJed Brown PetscReal mu; 25c4762a1bSJed Brown PetscReal next_output; 26c4762a1bSJed Brown }; 27c4762a1bSJed Brown 28c4762a1bSJed Brown /* 290e3d61c9SBarry Smith User-defined routines 30c4762a1bSJed Brown */ 31*2a8381b2SBarry Smith static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec X, Vec F, PetscCtx ctx) 32d71ae5a4SJacob Faibussowitsch { 33c4762a1bSJed Brown User user = (User)ctx; 34c4762a1bSJed Brown PetscScalar *f; 35c4762a1bSJed Brown const PetscScalar *x; 36c4762a1bSJed Brown 37c4762a1bSJed Brown PetscFunctionBeginUser; 389566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 399566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 40c4762a1bSJed Brown f[0] = x[1]; 41c4762a1bSJed Brown f[1] = user->mu * (1. - x[0] * x[0]) * x[1] - x[0]; 429566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 439566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 45c4762a1bSJed Brown } 46c4762a1bSJed Brown 47*2a8381b2SBarry Smith static PetscErrorCode IFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, PetscCtx ctx) 48d71ae5a4SJacob Faibussowitsch { 49c4762a1bSJed Brown User user = (User)ctx; 50c4762a1bSJed Brown const PetscScalar *x, *xdot; 51c4762a1bSJed Brown PetscScalar *f; 52c4762a1bSJed Brown 53c4762a1bSJed Brown PetscFunctionBeginUser; 549566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 559566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Xdot, &xdot)); 569566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 57c4762a1bSJed Brown f[0] = xdot[0] - x[1]; 58c4762a1bSJed Brown f[1] = xdot[1] - user->mu * ((1.0 - x[0] * x[0]) * x[1] - x[0]); 599566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 609566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Xdot, &xdot)); 619566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 63c4762a1bSJed Brown } 64c4762a1bSJed Brown 65*2a8381b2SBarry Smith static PetscErrorCode IJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat A, Mat B, PetscCtx ctx) 66d71ae5a4SJacob Faibussowitsch { 67c4762a1bSJed Brown User user = (User)ctx; 68c4762a1bSJed Brown PetscInt rowcol[] = {0, 1}; 69c4762a1bSJed Brown const PetscScalar *x; 70c4762a1bSJed Brown PetscScalar J[2][2]; 71c4762a1bSJed Brown 72c4762a1bSJed Brown PetscFunctionBeginUser; 739566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 749371c9d4SSatish Balay J[0][0] = a; 759371c9d4SSatish Balay J[0][1] = -1.0; 769371c9d4SSatish Balay J[1][0] = user->mu * (2.0 * x[0] * x[1] + 1.0); 779371c9d4SSatish Balay J[1][1] = a - user->mu * (1.0 - x[0] * x[0]); 789566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES)); 799566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 80c4762a1bSJed Brown 819566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 829566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 83c4762a1bSJed Brown if (A != B) { 849566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 859566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 86c4762a1bSJed Brown } 873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 88c4762a1bSJed Brown } 89c4762a1bSJed Brown 90c4762a1bSJed Brown /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */ 91*2a8381b2SBarry Smith static PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal t, Vec X, PetscCtx ctx) 92d71ae5a4SJacob Faibussowitsch { 93c4762a1bSJed Brown const PetscScalar *x; 94c4762a1bSJed Brown PetscReal tfinal, dt; 95c4762a1bSJed Brown User user = (User)ctx; 96c4762a1bSJed Brown Vec interpolatedX; 97c4762a1bSJed Brown 98c4762a1bSJed Brown PetscFunctionBeginUser; 999566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts, &dt)); 1009566063dSJacob Faibussowitsch PetscCall(TSGetMaxTime(ts, &tfinal)); 101c4762a1bSJed Brown 102c4762a1bSJed Brown while (user->next_output <= t && user->next_output <= tfinal) { 1039566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X, &interpolatedX)); 1049566063dSJacob Faibussowitsch PetscCall(TSInterpolate(ts, user->next_output, interpolatedX)); 1059566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(interpolatedX, &x)); 1069371c9d4SSatish Balay 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]))); 1079566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(interpolatedX, &x)); 1089566063dSJacob Faibussowitsch PetscCall(VecDestroy(&interpolatedX)); 109c4762a1bSJed Brown user->next_output += 0.1; 110c4762a1bSJed Brown } 1113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 112c4762a1bSJed Brown } 113c4762a1bSJed Brown 114d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 115d71ae5a4SJacob Faibussowitsch { 116c4762a1bSJed Brown TS ts; /* nonlinear solver */ 117c4762a1bSJed Brown Vec x; /* solution, residual vectors */ 118c4762a1bSJed Brown Mat A; /* Jacobian matrix */ 119c4762a1bSJed Brown PetscInt steps; 120c4762a1bSJed Brown PetscReal ftime = 0.5; 121c4762a1bSJed Brown PetscBool monitor = PETSC_FALSE, implicitform = PETSC_TRUE; 122c4762a1bSJed Brown PetscScalar *x_ptr; 123c4762a1bSJed Brown PetscMPIInt size; 124c4762a1bSJed Brown struct _n_User user; 125c4762a1bSJed Brown 126c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 127c4762a1bSJed Brown Initialize program 128c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 129327415f7SBarry Smith PetscFunctionBeginUser; 1309566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 1319566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 1323c633725SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!"); 133c4762a1bSJed Brown 134c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 135c4762a1bSJed Brown Set runtime options 136c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 137c4762a1bSJed Brown user.next_output = 0.0; 138c4762a1bSJed Brown user.mu = 1.0e3; 1399566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-monitor", &monitor, NULL)); 1409566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-implicitform", &implicitform, NULL)); 141d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Physical parameters", NULL); 1429566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-mu", "Stiffness parameter", "<1.0e6>", user.mu, &user.mu, NULL)); 143d0609cedSBarry Smith PetscOptionsEnd(); 144c4762a1bSJed Brown 145c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 146c4762a1bSJed Brown Create necessary matrix and vectors, solve same ODE on every process 147c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1489566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 1499566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, 2, 2)); 1509566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 1519566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 152c4762a1bSJed Brown 1539566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &x, NULL)); 154c4762a1bSJed Brown 155c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 156c4762a1bSJed Brown Create timestepping solver context 157c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1589566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 159c4762a1bSJed Brown if (implicitform) { 1609566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts, NULL, IFunction, &user)); 1619566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts, A, A, IJacobian, &user)); 1629566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSBEULER)); 163c4762a1bSJed Brown } else { 1649566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, &user)); 1659566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSRK)); 166c4762a1bSJed Brown } 1679566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, ftime)); 1689566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, 0.001)); 1699566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER)); 17048a46eb9SPierre Jolivet if (monitor) PetscCall(TSMonitorSet(ts, Monitor, &user, NULL)); 171c4762a1bSJed Brown 172c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 173c4762a1bSJed Brown Set initial conditions 174c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1759566063dSJacob Faibussowitsch PetscCall(VecGetArray(x, &x_ptr)); 176c4762a1bSJed Brown x_ptr[0] = 2.0; 177c4762a1bSJed Brown x_ptr[1] = -2.0 / 3.0 + 10.0 / (81.0 * user.mu) - 292.0 / (2187.0 * user.mu * user.mu); 1789566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x, &x_ptr)); 179c4762a1bSJed Brown 180c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 181c4762a1bSJed Brown Set runtime options 182c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1839566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 184c4762a1bSJed Brown 185c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 186c4762a1bSJed Brown Solve nonlinear system 187c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1889566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, x)); 1899566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts, &ftime)); 1909566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &steps)); 19163a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "steps %" PetscInt_FMT ", ftime %g\n", steps, (double)ftime)); 1929566063dSJacob Faibussowitsch PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD)); 193c4762a1bSJed Brown 194c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 195c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 196c4762a1bSJed Brown are no longer needed. 197c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1989566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1999566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 2009566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 201c4762a1bSJed Brown 2029566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 2034ad8454bSPierre Jolivet return 0; 204c4762a1bSJed Brown } 205c4762a1bSJed Brown 206c4762a1bSJed Brown /*TEST 207c4762a1bSJed Brown 208c4762a1bSJed Brown test: 209c4762a1bSJed Brown requires: !single 210c4762a1bSJed Brown args: -mu 1e6 211c4762a1bSJed Brown 21210b8587eSHendrik Ranocha test: 21310b8587eSHendrik Ranocha requires: !single 21410b8587eSHendrik Ranocha suffix: 2 21510b8587eSHendrik Ranocha args: -implicitform false -ts_type rk -ts_rk_type 5dp -ts_adapt_type dsp 21610b8587eSHendrik Ranocha 21710b8587eSHendrik Ranocha test: 21810b8587eSHendrik Ranocha requires: !single 21910b8587eSHendrik Ranocha suffix: 3 22010b8587eSHendrik Ranocha args: -implicitform false -ts_type rk -ts_rk_type 5dp -ts_adapt_type dsp -ts_adapt_dsp_filter H0312 22110b8587eSHendrik Ranocha 222c4762a1bSJed Brown TEST*/ 223