1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Solves the van der Pol equation.\n\ 3c4762a1bSJed Brown Input parameters include:\n"; 4c4762a1bSJed Brown 5c4762a1bSJed Brown /* ------------------------------------------------------------------------ 6c4762a1bSJed Brown 7c4762a1bSJed Brown This program solves the van der Pol DAE ODE equivalent 8c4762a1bSJed Brown y' = z (1) 9c4762a1bSJed Brown z' = \mu ((1-y^2)z-y) 10c4762a1bSJed Brown on the domain 0 <= x <= 1, with the boundary conditions 11c4762a1bSJed Brown y(0) = 2, y'(0) = - 2/3 +10/(81*\mu) - 292/(2187*\mu^2), 12c4762a1bSJed Brown and 13c4762a1bSJed Brown \mu = 10^6 ( y'(0) ~ -0.6666665432100101). 14c4762a1bSJed 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. 15c4762a1bSJed Brown 16c4762a1bSJed Brown Notes: 17c4762a1bSJed Brown This code demonstrates the TS solver interface to an ODE -- RHSFunction for explicit form and IFunction for implicit form. 18c4762a1bSJed Brown 19c4762a1bSJed Brown ------------------------------------------------------------------------- */ 20c4762a1bSJed Brown 21c4762a1bSJed Brown #include <petscts.h> 22c4762a1bSJed Brown 23c4762a1bSJed Brown typedef struct _n_User *User; 24c4762a1bSJed Brown struct _n_User { 25c4762a1bSJed Brown PetscReal mu; 26c4762a1bSJed Brown PetscReal next_output; 27c4762a1bSJed Brown }; 28c4762a1bSJed Brown 29c4762a1bSJed Brown /* 300e3d61c9SBarry Smith User-defined routines 31c4762a1bSJed Brown */ 32d71ae5a4SJacob Faibussowitsch static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec X, Vec F, void *ctx) 33d71ae5a4SJacob Faibussowitsch { 34c4762a1bSJed Brown User user = (User)ctx; 35c4762a1bSJed Brown PetscScalar *f; 36c4762a1bSJed Brown const PetscScalar *x; 37c4762a1bSJed Brown 38c4762a1bSJed Brown PetscFunctionBeginUser; 399566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 409566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 41c4762a1bSJed Brown f[0] = x[1]; 42c4762a1bSJed Brown f[1] = user->mu * (1. - x[0] * x[0]) * x[1] - x[0]; 439566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 449566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 45*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 46c4762a1bSJed Brown } 47c4762a1bSJed Brown 48d71ae5a4SJacob Faibussowitsch static PetscErrorCode IFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, void *ctx) 49d71ae5a4SJacob Faibussowitsch { 50c4762a1bSJed Brown User user = (User)ctx; 51c4762a1bSJed Brown const PetscScalar *x, *xdot; 52c4762a1bSJed Brown PetscScalar *f; 53c4762a1bSJed Brown 54c4762a1bSJed Brown PetscFunctionBeginUser; 559566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 569566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Xdot, &xdot)); 579566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 58c4762a1bSJed Brown f[0] = xdot[0] - x[1]; 59c4762a1bSJed Brown f[1] = xdot[1] - user->mu * ((1.0 - x[0] * x[0]) * x[1] - x[0]); 609566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 619566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Xdot, &xdot)); 629566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 63*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 64c4762a1bSJed Brown } 65c4762a1bSJed Brown 66d71ae5a4SJacob Faibussowitsch static PetscErrorCode IJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat A, Mat B, void *ctx) 67d71ae5a4SJacob Faibussowitsch { 68c4762a1bSJed Brown User user = (User)ctx; 69c4762a1bSJed Brown PetscInt rowcol[] = {0, 1}; 70c4762a1bSJed Brown const PetscScalar *x; 71c4762a1bSJed Brown PetscScalar J[2][2]; 72c4762a1bSJed Brown 73c4762a1bSJed Brown PetscFunctionBeginUser; 749566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 759371c9d4SSatish Balay J[0][0] = a; 769371c9d4SSatish Balay J[0][1] = -1.0; 779371c9d4SSatish Balay J[1][0] = user->mu * (2.0 * x[0] * x[1] + 1.0); 789371c9d4SSatish Balay J[1][1] = a - user->mu * (1.0 - x[0] * x[0]); 799566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES)); 809566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 81c4762a1bSJed Brown 829566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 839566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 84c4762a1bSJed Brown if (A != B) { 859566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 869566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 87c4762a1bSJed Brown } 88*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 89c4762a1bSJed Brown } 90c4762a1bSJed Brown 91c4762a1bSJed Brown /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */ 92d71ae5a4SJacob Faibussowitsch static PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal t, Vec X, void *ctx) 93d71ae5a4SJacob Faibussowitsch { 94c4762a1bSJed Brown const PetscScalar *x; 95c4762a1bSJed Brown PetscReal tfinal, dt; 96c4762a1bSJed Brown User user = (User)ctx; 97c4762a1bSJed Brown Vec interpolatedX; 98c4762a1bSJed Brown 99c4762a1bSJed Brown PetscFunctionBeginUser; 1009566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts, &dt)); 1019566063dSJacob Faibussowitsch PetscCall(TSGetMaxTime(ts, &tfinal)); 102c4762a1bSJed Brown 103c4762a1bSJed Brown while (user->next_output <= t && user->next_output <= tfinal) { 1049566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X, &interpolatedX)); 1059566063dSJacob Faibussowitsch PetscCall(TSInterpolate(ts, user->next_output, interpolatedX)); 1069566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(interpolatedX, &x)); 1079371c9d4SSatish 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]))); 1089566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(interpolatedX, &x)); 1099566063dSJacob Faibussowitsch PetscCall(VecDestroy(&interpolatedX)); 110c4762a1bSJed Brown user->next_output += 0.1; 111c4762a1bSJed Brown } 112*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 113c4762a1bSJed Brown } 114c4762a1bSJed Brown 115d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 116d71ae5a4SJacob Faibussowitsch { 117c4762a1bSJed Brown TS ts; /* nonlinear solver */ 118c4762a1bSJed Brown Vec x; /* solution, residual vectors */ 119c4762a1bSJed Brown Mat A; /* Jacobian matrix */ 120c4762a1bSJed Brown PetscInt steps; 121c4762a1bSJed Brown PetscReal ftime = 0.5; 122c4762a1bSJed Brown PetscBool monitor = PETSC_FALSE, implicitform = PETSC_TRUE; 123c4762a1bSJed Brown PetscScalar *x_ptr; 124c4762a1bSJed Brown PetscMPIInt size; 125c4762a1bSJed Brown struct _n_User user; 126c4762a1bSJed Brown 127c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 128c4762a1bSJed Brown Initialize program 129c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 130327415f7SBarry Smith PetscFunctionBeginUser; 1319566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 1329566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 1333c633725SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!"); 134c4762a1bSJed Brown 135c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 136c4762a1bSJed Brown Set runtime options 137c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 138c4762a1bSJed Brown user.next_output = 0.0; 139c4762a1bSJed Brown user.mu = 1.0e3; 1409566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-monitor", &monitor, NULL)); 1419566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-implicitform", &implicitform, NULL)); 142d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Physical parameters", NULL); 1439566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-mu", "Stiffness parameter", "<1.0e6>", user.mu, &user.mu, NULL)); 144d0609cedSBarry Smith PetscOptionsEnd(); 145c4762a1bSJed Brown 146c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 147c4762a1bSJed Brown Create necessary matrix and vectors, solve same ODE on every process 148c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1499566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 1509566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, 2, 2)); 1519566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 1529566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 153c4762a1bSJed Brown 1549566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &x, NULL)); 155c4762a1bSJed Brown 156c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 157c4762a1bSJed Brown Create timestepping solver context 158c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1599566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 160c4762a1bSJed Brown if (implicitform) { 1619566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts, NULL, IFunction, &user)); 1629566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts, A, A, IJacobian, &user)); 1639566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSBEULER)); 164c4762a1bSJed Brown } else { 1659566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, &user)); 1669566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSRK)); 167c4762a1bSJed Brown } 1689566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, ftime)); 1699566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, 0.001)); 1709566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER)); 17148a46eb9SPierre Jolivet if (monitor) PetscCall(TSMonitorSet(ts, Monitor, &user, NULL)); 172c4762a1bSJed Brown 173c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 174c4762a1bSJed Brown Set initial conditions 175c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1769566063dSJacob Faibussowitsch PetscCall(VecGetArray(x, &x_ptr)); 177c4762a1bSJed Brown x_ptr[0] = 2.0; 178c4762a1bSJed Brown x_ptr[1] = -2.0 / 3.0 + 10.0 / (81.0 * user.mu) - 292.0 / (2187.0 * user.mu * user.mu); 1799566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x, &x_ptr)); 180c4762a1bSJed Brown 181c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 182c4762a1bSJed Brown Set runtime options 183c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1849566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 185c4762a1bSJed Brown 186c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 187c4762a1bSJed Brown Solve nonlinear system 188c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1899566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, x)); 1909566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts, &ftime)); 1919566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &steps)); 19263a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "steps %" PetscInt_FMT ", ftime %g\n", steps, (double)ftime)); 1939566063dSJacob Faibussowitsch PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD)); 194c4762a1bSJed Brown 195c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 196c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 197c4762a1bSJed Brown are no longer needed. 198c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1999566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 2009566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 2019566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 202c4762a1bSJed Brown 2039566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 204d0609cedSBarry Smith return (0); 205c4762a1bSJed Brown } 206c4762a1bSJed Brown 207c4762a1bSJed Brown /*TEST 208c4762a1bSJed Brown 209c4762a1bSJed Brown test: 210c4762a1bSJed Brown requires: !single 211c4762a1bSJed Brown args: -mu 1e6 212c4762a1bSJed Brown 21310b8587eSHendrik Ranocha test: 21410b8587eSHendrik Ranocha requires: !single 21510b8587eSHendrik Ranocha suffix: 2 21610b8587eSHendrik Ranocha args: -implicitform false -ts_type rk -ts_rk_type 5dp -ts_adapt_type dsp 21710b8587eSHendrik Ranocha 21810b8587eSHendrik Ranocha test: 21910b8587eSHendrik Ranocha requires: !single 22010b8587eSHendrik Ranocha suffix: 3 22110b8587eSHendrik Ranocha args: -implicitform false -ts_type rk -ts_rk_type 5dp -ts_adapt_type dsp -ts_adapt_dsp_filter H0312 22210b8587eSHendrik Ranocha 223c4762a1bSJed Brown TEST*/ 224