1c4762a1bSJed Brown static char help[] = "Solves the van der Pol DAE.\n\ 2c4762a1bSJed Brown Input parameters include:\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown /* ------------------------------------------------------------------------ 5c4762a1bSJed Brown 6c4762a1bSJed Brown This program solves the van der Pol DAE 7c4762a1bSJed Brown y' = -z = f(y,z) (1) 8c4762a1bSJed Brown 0 = y-(z^3/3 - z) = g(y,z) 9c4762a1bSJed Brown on the domain 0 <= x <= 1, with the boundary conditions 10c4762a1bSJed Brown y(0) = -2, y'(0) = -2.355301397608119909925287735864250951918 11c4762a1bSJed Brown This is a nonlinear equation. 12c4762a1bSJed Brown 13c4762a1bSJed Brown Notes: 14c4762a1bSJed Brown This code demonstrates the TS solver interface with the Van der Pol DAE, 15c4762a1bSJed Brown namely it is the case when there is no RHS (meaning the RHS == 0), and the 16c4762a1bSJed Brown equations are converted to two variants of linear problems, u_t = f(u,t), 17c4762a1bSJed Brown namely turning (1) into a vector equation in terms of u, 18c4762a1bSJed Brown 19c4762a1bSJed Brown [ y' + z ] = [ 0 ] 20c4762a1bSJed Brown [ (z^3/3 - z) - y ] [ 0 ] 21c4762a1bSJed Brown 22c4762a1bSJed Brown which then we can write as a vector equation 23c4762a1bSJed Brown 24c4762a1bSJed Brown [ u_1' + u_2 ] = [ 0 ] (2) 25c4762a1bSJed Brown [ (u_2^3/3 - u_2) - u_1 ] [ 0 ] 26c4762a1bSJed Brown 27c4762a1bSJed Brown which is now in the desired form of u_t = f(u,t). As this is a DAE, and 28c4762a1bSJed Brown there is no u_2', there is no need for a split, 29c4762a1bSJed Brown 30c4762a1bSJed Brown so 31c4762a1bSJed Brown 325ab1ac2bSHong Zhang [ F(u',u,t) ] = [ u_1' ] + [ u_2 ] 33c4762a1bSJed Brown [ 0 ] [ (u_2^3/3 - u_2) - u_1 ] 34c4762a1bSJed Brown 355ab1ac2bSHong Zhang Using the definition of the Jacobian of F (from the PETSc user manual), 365ab1ac2bSHong Zhang in the equation F(u',u,t) = G(u,t), 37c4762a1bSJed Brown 385ab1ac2bSHong Zhang dF dF 395ab1ac2bSHong Zhang J(F) = a * -- - -- 40c4762a1bSJed Brown du' du 41c4762a1bSJed Brown 42c4762a1bSJed Brown where d is the partial derivative. In this example, 43c4762a1bSJed Brown 445ab1ac2bSHong Zhang dF [ 1 ; 0 ] 45c4762a1bSJed Brown -- = [ ] 46c4762a1bSJed Brown du' [ 0 ; 0 ] 47c4762a1bSJed Brown 485ab1ac2bSHong Zhang dF [ 0 ; 1 ] 49c4762a1bSJed Brown -- = [ ] 50c4762a1bSJed Brown du [ -1 ; 1 - u_2^2 ] 51c4762a1bSJed Brown 52c4762a1bSJed Brown Hence, 53c4762a1bSJed Brown 54c4762a1bSJed Brown [ a ; -1 ] 555ab1ac2bSHong Zhang J(F) = [ ] 56c4762a1bSJed Brown [ 1 ; u_2^2 - 1 ] 57c4762a1bSJed Brown 58c4762a1bSJed Brown ------------------------------------------------------------------------- */ 59c4762a1bSJed Brown 60c4762a1bSJed Brown #include <petscts.h> 61c4762a1bSJed Brown 62c4762a1bSJed Brown typedef struct _n_User *User; 63c4762a1bSJed Brown struct _n_User { 64c4762a1bSJed Brown PetscReal next_output; 65c4762a1bSJed Brown }; 66c4762a1bSJed Brown 67c4762a1bSJed Brown /* 680e3d61c9SBarry Smith User-defined routines 69c4762a1bSJed Brown */ 70c4762a1bSJed Brown 71*2a8381b2SBarry Smith static PetscErrorCode IFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, PetscCtx ctx) 72d71ae5a4SJacob Faibussowitsch { 73c4762a1bSJed Brown PetscScalar *f; 74c4762a1bSJed Brown const PetscScalar *x, *xdot; 75c4762a1bSJed Brown 76c4762a1bSJed Brown PetscFunctionBeginUser; 779566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 789566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Xdot, &xdot)); 799566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 80c4762a1bSJed Brown f[0] = xdot[0] + x[1]; 81c4762a1bSJed Brown f[1] = (x[1] * x[1] * x[1] / 3.0 - x[1]) - x[0]; 829566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 839566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Xdot, &xdot)); 849566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 86c4762a1bSJed Brown } 87c4762a1bSJed Brown 88*2a8381b2SBarry Smith static PetscErrorCode IJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat A, Mat B, PetscCtx ctx) 89d71ae5a4SJacob Faibussowitsch { 90c4762a1bSJed Brown PetscInt rowcol[] = {0, 1}; 91c4762a1bSJed Brown PetscScalar J[2][2]; 92c4762a1bSJed Brown const PetscScalar *x; 93c4762a1bSJed Brown 94c4762a1bSJed Brown PetscFunctionBeginUser; 959566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 969371c9d4SSatish Balay J[0][0] = a; 979371c9d4SSatish Balay J[0][1] = -1.; 989371c9d4SSatish Balay J[1][0] = 1.; 999371c9d4SSatish Balay J[1][1] = -1. + x[1] * x[1]; 1009566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES)); 1019566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 102c4762a1bSJed Brown 1039566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 1049566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 105c4762a1bSJed Brown if (A != B) { 1069566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 1079566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 108c4762a1bSJed Brown } 1093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 110c4762a1bSJed Brown } 111c4762a1bSJed Brown 112d71ae5a4SJacob Faibussowitsch static PetscErrorCode RegisterMyARK2(void) 113d71ae5a4SJacob Faibussowitsch { 114c4762a1bSJed Brown PetscFunctionBeginUser; 115c4762a1bSJed Brown { 1169371c9d4SSatish Balay const PetscReal A[3][3] = 1179371c9d4SSatish Balay { 1189371c9d4SSatish Balay {0, 0, 0}, 119c4762a1bSJed Brown {0.41421356237309504880, 0, 0}, 1209371c9d4SSatish Balay {0.75, 0.25, 0} 1219371c9d4SSatish Balay }, 1229371c9d4SSatish Balay At[3][3] = {{0, 0, 0}, {0.12132034355964257320, 0.29289321881345247560, 0}, {0.20710678118654752440, 0.50000000000000000000, 0.29289321881345247560}}, *bembedt = NULL, *bembed = NULL; 1239566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister("myark2", 2, 3, &At[0][0], NULL, NULL, &A[0][0], NULL, NULL, bembedt, bembed, 0, NULL, NULL)); 124c4762a1bSJed Brown } 1253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 126c4762a1bSJed Brown } 127c4762a1bSJed Brown 128c4762a1bSJed Brown /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */ 129*2a8381b2SBarry Smith static PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal t, Vec X, PetscCtx ctx) 130d71ae5a4SJacob Faibussowitsch { 131c4762a1bSJed Brown const PetscScalar *x; 132c4762a1bSJed Brown PetscReal tfinal, dt; 133c4762a1bSJed Brown User user = (User)ctx; 134c4762a1bSJed Brown Vec interpolatedX; 135c4762a1bSJed Brown 136c4762a1bSJed Brown PetscFunctionBeginUser; 1379566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts, &dt)); 1389566063dSJacob Faibussowitsch PetscCall(TSGetMaxTime(ts, &tfinal)); 139c4762a1bSJed Brown 140c4762a1bSJed Brown while (user->next_output <= t && user->next_output <= tfinal) { 1419566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X, &interpolatedX)); 1429566063dSJacob Faibussowitsch PetscCall(TSInterpolate(ts, user->next_output, interpolatedX)); 1439566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(interpolatedX, &x)); 14463a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "[%.1f] %3" 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]))); 1459566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(interpolatedX, &x)); 1469566063dSJacob Faibussowitsch PetscCall(VecDestroy(&interpolatedX)); 147c4762a1bSJed Brown user->next_output += PetscRealConstant(0.1); 148c4762a1bSJed Brown } 1493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 150c4762a1bSJed Brown } 151c4762a1bSJed Brown 152d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 153d71ae5a4SJacob Faibussowitsch { 154c4762a1bSJed Brown TS ts; /* nonlinear solver */ 155c4762a1bSJed Brown Vec x; /* solution, residual vectors */ 156c4762a1bSJed Brown Mat A; /* Jacobian matrix */ 157c4762a1bSJed Brown PetscInt steps; 158c4762a1bSJed Brown PetscReal ftime = 0.5; 159c4762a1bSJed Brown PetscBool monitor = PETSC_FALSE; 160c4762a1bSJed Brown PetscScalar *x_ptr; 161c4762a1bSJed Brown PetscMPIInt size; 162c4762a1bSJed Brown struct _n_User user; 163c4762a1bSJed Brown 164c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 165c4762a1bSJed Brown Initialize program 166c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 167327415f7SBarry Smith PetscFunctionBeginUser; 1689566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 1699566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 1703c633725SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!"); 171c4762a1bSJed Brown 1729566063dSJacob Faibussowitsch PetscCall(RegisterMyARK2()); 173c4762a1bSJed Brown 174c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 175c4762a1bSJed Brown Set runtime options 176c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 177c4762a1bSJed Brown 178c4762a1bSJed Brown user.next_output = 0.0; 1799566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-monitor", &monitor, NULL)); 180c4762a1bSJed Brown 181c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 182c4762a1bSJed Brown Create necessary matrix and vectors, solve same ODE on every process 183c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1849566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 1859566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, 2, 2)); 1869566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 1879566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 1889566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &x, NULL)); 189c4762a1bSJed Brown 190c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 191c4762a1bSJed Brown Create timestepping solver context 192c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1939566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 1949566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSBEULER)); 1959566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts, NULL, IFunction, &user)); 1969566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts, A, A, IJacobian, &user)); 1979566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, ftime)); 1989566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER)); 19948a46eb9SPierre Jolivet if (monitor) PetscCall(TSMonitorSet(ts, Monitor, &user, NULL)); 200c4762a1bSJed Brown 201c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 202c4762a1bSJed Brown Set initial conditions 203c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2049566063dSJacob Faibussowitsch PetscCall(VecGetArray(x, &x_ptr)); 2059371c9d4SSatish Balay x_ptr[0] = -2; 2069371c9d4SSatish Balay x_ptr[1] = -2.355301397608119909925287735864250951918; 2079566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x, &x_ptr)); 2089566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, .001)); 209c4762a1bSJed Brown 210c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 211c4762a1bSJed Brown Set runtime options 212c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2139566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 214c4762a1bSJed Brown 215c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 216c4762a1bSJed Brown Solve nonlinear system 217c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2189566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, x)); 2199566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts, &ftime)); 2209566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &steps)); 22163a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "steps %3" PetscInt_FMT ", ftime %g\n", steps, (double)ftime)); 2229566063dSJacob Faibussowitsch PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD)); 223c4762a1bSJed Brown 224c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 225c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 226c4762a1bSJed Brown are no longer needed. 227c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2289566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 2299566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 2309566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 231c4762a1bSJed Brown 2329566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 233b122ec5aSJacob Faibussowitsch return 0; 234c4762a1bSJed Brown } 235c4762a1bSJed Brown 236c4762a1bSJed Brown /*TEST 237c4762a1bSJed Brown 238c4762a1bSJed Brown test: 239c4762a1bSJed Brown requires: !single 240c4762a1bSJed Brown suffix: a 241c4762a1bSJed Brown args: -monitor -ts_type bdf -ts_rtol 1e-6 -ts_atol 1e-6 -ts_view -ts_adapt_type dsp 242c4762a1bSJed Brown output_file: output/ex19_pi42.out 243c4762a1bSJed Brown 244c4762a1bSJed Brown test: 245c4762a1bSJed Brown requires: !single 246c4762a1bSJed Brown suffix: b 247c4762a1bSJed Brown args: -monitor -ts_type bdf -ts_rtol 1e-6 -ts_atol 1e-6 -ts_view -ts_adapt_type dsp -ts_adapt_dsp_filter PI42 248c4762a1bSJed Brown output_file: output/ex19_pi42.out 249c4762a1bSJed Brown 250c4762a1bSJed Brown test: 251c4762a1bSJed Brown requires: !single 252c4762a1bSJed Brown suffix: c 253c4762a1bSJed Brown args: -monitor -ts_type bdf -ts_rtol 1e-6 -ts_atol 1e-6 -ts_view -ts_adapt_type dsp -ts_adapt_dsp_pid 0.4,0.2 254c4762a1bSJed Brown output_file: output/ex19_pi42.out 255c4762a1bSJed Brown 256e5b8ffdfSLisandro Dalcin test: 257e5b8ffdfSLisandro Dalcin requires: !single 258e5b8ffdfSLisandro Dalcin suffix: bdf_reject 259188af4bfSBarry Smith args: -ts_type bdf -ts_time_step 0.5 -ts_max_steps 1 -ts_max_step_rejections {{0 1 2}separate_output} -ts_error_if_step_fails false -ts_adapt_monitor 260e5b8ffdfSLisandro Dalcin 261c4762a1bSJed Brown TEST*/ 262