1c4762a1bSJed Brown /* 2c4762a1bSJed Brown * ex_vdp.c 3c4762a1bSJed Brown * 4c4762a1bSJed Brown * Created on: Jun 1, 2012 5c4762a1bSJed Brown * Author: Hong Zhang 6c4762a1bSJed Brown */ 7c4762a1bSJed Brown static char help[] = "Solves the van der Pol equation. \n Input parameters include:\n"; 8c4762a1bSJed Brown 9c4762a1bSJed Brown /* 10c4762a1bSJed Brown * This program solves the van der Pol equation 11c4762a1bSJed Brown * y' = z (1) 12c4762a1bSJed Brown * z' = (((1-y^2)*z-y)/eps (2) 13c4762a1bSJed Brown * on the domain 0<=x<=0.5, with the initial conditions 14c4762a1bSJed Brown * y(0) = 2, 15c4762a1bSJed Brown * z(0) = -2/3 + 10/81*eps - 292/2187*eps^2-1814/19683*eps^3 16c4762a1bSJed Brown * IMEX schemes are applied to the splitted equation 17c4762a1bSJed Brown * [y'] = [z] + [0 ] 18c4762a1bSJed Brown * [z'] [0] [(((1-y^2)*z-y)/eps] 19c4762a1bSJed Brown * 20c4762a1bSJed Brown * F(x)= [z] 21c4762a1bSJed Brown * [0] 22c4762a1bSJed Brown * 23c4762a1bSJed Brown * G(x)= [y'] - [0 ] 24c4762a1bSJed Brown * [z'] [(((1-y^2)*z-y)/eps] 25c4762a1bSJed Brown * 26c4762a1bSJed Brown * JG(x) = G_x + a G_xdot 27c4762a1bSJed Brown */ 28c4762a1bSJed Brown 29c4762a1bSJed Brown #include <petscdmda.h> 30c4762a1bSJed Brown #include <petscts.h> 31c4762a1bSJed Brown 32c4762a1bSJed Brown typedef struct _User *User; 33c4762a1bSJed Brown struct _User { 34c4762a1bSJed Brown PetscReal mu; /*stiffness control coefficient: epsilon*/ 35c4762a1bSJed Brown }; 36c4762a1bSJed Brown 37c4762a1bSJed Brown static PetscErrorCode RHSFunction(TS, PetscReal, Vec, Vec, void *); 38c4762a1bSJed Brown static PetscErrorCode IFunction(TS, PetscReal, Vec, Vec, Vec, void *); 39c4762a1bSJed Brown static PetscErrorCode IJacobian(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *); 40c4762a1bSJed Brown 41*9371c9d4SSatish Balay int main(int argc, char **argv) { 42c4762a1bSJed Brown TS ts; 43c4762a1bSJed Brown Vec x; /* solution vector */ 44c4762a1bSJed Brown Mat A; /* Jacobian */ 45c4762a1bSJed Brown PetscInt steps, mx, eimex_rowcol[2], two; 46c4762a1bSJed Brown PetscScalar *x_ptr; 47c4762a1bSJed Brown PetscReal ftime, dt, norm; 48c4762a1bSJed Brown Vec ref; 49c4762a1bSJed Brown struct _User user; /* user-defined work context */ 50c4762a1bSJed Brown PetscViewer viewer; 51c4762a1bSJed Brown 52327415f7SBarry Smith PetscFunctionBeginUser; 539566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 54c4762a1bSJed Brown /* Initialize user application context */ 55d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "van der Pol options", ""); 56c4762a1bSJed Brown user.mu = 1e0; 579566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-eps", "Stiffness controller", "", user.mu, &user.mu, NULL)); 58d0609cedSBarry Smith PetscOptionsEnd(); 59c4762a1bSJed Brown 60c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 61c4762a1bSJed Brown Set runtime options 62c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 63c4762a1bSJed Brown /* 649566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL)); 65c4762a1bSJed Brown */ 66c4762a1bSJed Brown 67c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 68c4762a1bSJed Brown Create necessary matrix and vectors, solve same ODE on every process 69c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 709566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 719566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, 2, 2)); 729566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 739566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 749566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &x, NULL)); 75c4762a1bSJed Brown 769566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &ref, NULL)); 779566063dSJacob Faibussowitsch PetscCall(VecGetArray(ref, &x_ptr)); 78c4762a1bSJed Brown /* 79c4762a1bSJed Brown * [0,1], mu=10^-3 80c4762a1bSJed Brown */ 81c4762a1bSJed Brown /* 82c4762a1bSJed Brown x_ptr[0] = -1.8881254106283; 83c4762a1bSJed Brown x_ptr[1] = 0.7359074233370;*/ 84c4762a1bSJed Brown 85c4762a1bSJed Brown /* 86c4762a1bSJed Brown * [0,0.5],mu=10^-3 87c4762a1bSJed Brown */ 88c4762a1bSJed Brown /* 89c4762a1bSJed Brown x_ptr[0] = 1.596980778659137; 90c4762a1bSJed Brown x_ptr[1] = -1.029103015879544; 91c4762a1bSJed Brown */ 92c4762a1bSJed Brown /* 93c4762a1bSJed Brown * [0,0.5],mu=1 94c4762a1bSJed Brown */ 95c4762a1bSJed Brown x_ptr[0] = 1.619084329683235; 96c4762a1bSJed Brown x_ptr[1] = -0.803530465176385; 97c4762a1bSJed Brown 98c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 99c4762a1bSJed Brown Create timestepping solver context 100c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1019566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 1029566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSEIMEX)); 1039566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, &user)); 1049566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts, NULL, IFunction, &user)); 1059566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts, A, A, IJacobian, &user)); 106c4762a1bSJed Brown 107c4762a1bSJed Brown dt = 0.00001; 108c4762a1bSJed Brown ftime = 1.1; 1099566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, dt)); 1109566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, ftime)); 1119566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER)); 112c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 113c4762a1bSJed Brown Set initial conditions 114c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1159566063dSJacob Faibussowitsch PetscCall(VecGetArray(x, &x_ptr)); 116c4762a1bSJed Brown x_ptr[0] = 2.; 117*9371c9d4SSatish Balay x_ptr[1] = -2. / 3. + 10. / 81. * (user.mu) - 292. / 2187. * (user.mu) * (user.mu) - 1814. / 19683. * (user.mu) * (user.mu) * (user.mu); 1189566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, x)); 1199566063dSJacob Faibussowitsch PetscCall(VecGetSize(x, &mx)); 120c4762a1bSJed Brown 121c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 122c4762a1bSJed Brown Set runtime options 123c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1249566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 125c4762a1bSJed Brown 126c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 127c4762a1bSJed Brown Solve nonlinear system 128c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1299566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, x)); 1309566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts, &ftime)); 1319566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &steps)); 132c4762a1bSJed Brown 1339566063dSJacob Faibussowitsch PetscCall(VecAXPY(x, -1.0, ref)); 1349566063dSJacob Faibussowitsch PetscCall(VecNorm(x, NORM_2, &norm)); 1359566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts, &dt)); 136c4762a1bSJed Brown 137*9371c9d4SSatish Balay eimex_rowcol[0] = 0; 138*9371c9d4SSatish Balay eimex_rowcol[1] = 0; 139*9371c9d4SSatish Balay two = 2; 1409566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-ts_eimex_row_col", eimex_rowcol, &two, NULL)); 1419566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "order %11s %18s %37s\n", "dt", "norm", "final solution components 0 and 1")); 1429566063dSJacob Faibussowitsch PetscCall(VecGetArray(x, &x_ptr)); 14363a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "(%" PetscInt_FMT ",%" PetscInt_FMT ") %10.8f %18.15f %18.15f %18.15f\n", eimex_rowcol[0], eimex_rowcol[1], (double)dt, (double)norm, (double)PetscRealPart(x_ptr[0]), (double)PetscRealPart(x_ptr[1]))); 1449566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x, &x_ptr)); 145c4762a1bSJed Brown 146c4762a1bSJed Brown /* Write line in convergence log */ 1479566063dSJacob Faibussowitsch PetscCall(PetscViewerCreate(PETSC_COMM_WORLD, &viewer)); 1489566063dSJacob Faibussowitsch PetscCall(PetscViewerSetType(viewer, PETSCVIEWERASCII)); 1499566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetMode(viewer, FILE_MODE_APPEND)); 1509566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetName(viewer, "eimex_nonstiff_vdp.txt")); 15163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " %" PetscInt_FMT " %10.8f %18.15f\n", eimex_rowcol[0], eimex_rowcol[1], (double)dt, (double)norm)); 1529566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 153c4762a1bSJed Brown 154c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 155c4762a1bSJed Brown Free work space. 156c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1579566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1589566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1599566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ref)); 1609566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 1619566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 162b122ec5aSJacob Faibussowitsch return 0; 163c4762a1bSJed Brown } 164c4762a1bSJed Brown 165*9371c9d4SSatish Balay static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec X, Vec F, void *ptr) { 166c4762a1bSJed Brown PetscScalar *f; 167c4762a1bSJed Brown const PetscScalar *x; 168c4762a1bSJed Brown 169c4762a1bSJed Brown PetscFunctionBegin; 1709566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 1719566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 172c4762a1bSJed Brown f[0] = x[1]; 173c4762a1bSJed Brown f[1] = 0.0; 1749566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 1759566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 176c4762a1bSJed Brown PetscFunctionReturn(0); 177c4762a1bSJed Brown } 178c4762a1bSJed Brown 179*9371c9d4SSatish Balay static PetscErrorCode IFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, void *ptr) { 180c4762a1bSJed Brown User user = (User)ptr; 181c4762a1bSJed Brown PetscScalar *f; 182c4762a1bSJed Brown const PetscScalar *x, *xdot; 183c4762a1bSJed Brown 184c4762a1bSJed Brown PetscFunctionBegin; 1859566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 1869566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Xdot, &xdot)); 1879566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 188c4762a1bSJed Brown f[0] = xdot[0]; 189c4762a1bSJed Brown f[1] = xdot[1] - ((1. - x[0] * x[0]) * x[1] - x[0]) / user->mu; 1909566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 1919566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Xdot, &xdot)); 1929566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 193c4762a1bSJed Brown PetscFunctionReturn(0); 194c4762a1bSJed Brown } 195c4762a1bSJed Brown 196*9371c9d4SSatish Balay static PetscErrorCode IJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat A, Mat B, void *ptr) { 197c4762a1bSJed Brown User user = (User)ptr; 198c4762a1bSJed Brown PetscReal mu = user->mu; 199c4762a1bSJed Brown PetscInt rowcol[] = {0, 1}; 200c4762a1bSJed Brown PetscScalar J[2][2]; 201c4762a1bSJed Brown const PetscScalar *x; 202c4762a1bSJed Brown 203c4762a1bSJed Brown PetscFunctionBegin; 2049566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 205c4762a1bSJed Brown J[0][0] = a; 206c4762a1bSJed Brown J[0][1] = 0; 207c4762a1bSJed Brown J[1][0] = (2. * x[0] * x[1] + 1.) / mu; 208c4762a1bSJed Brown J[1][1] = a - (1. - x[0] * x[0]) / mu; 2099566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES)); 2109566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 211c4762a1bSJed Brown 2129566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 2139566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 214c4762a1bSJed Brown if (A != B) { 2159566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 2169566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 217c4762a1bSJed Brown } 218c4762a1bSJed Brown PetscFunctionReturn(0); 219c4762a1bSJed Brown } 220c4762a1bSJed Brown 221c4762a1bSJed Brown /*TEST 222c4762a1bSJed Brown 223c4762a1bSJed Brown test: 224c4762a1bSJed Brown args: -ts_type eimex -ts_adapt_type none -pc_type lu -ts_dt 0.01 -ts_max_time 10 -ts_eimex_row_col 3,3 -ts_monitor_lg_solution 225c4762a1bSJed Brown requires: x 226c4762a1bSJed Brown 227c4762a1bSJed Brown test: 228c4762a1bSJed Brown suffix: adapt 229c4762a1bSJed Brown args: -ts_type eimex -ts_adapt_type none -pc_type lu -ts_dt 0.01 -ts_max_time 10 -ts_eimex_order_adapt -ts_eimex_max_rows 7 -ts_monitor_lg_solution 230c4762a1bSJed Brown requires: x 231c4762a1bSJed Brown 232c4762a1bSJed Brown test: 233c4762a1bSJed Brown suffix: loop 234c4762a1bSJed Brown args: -ts_type eimex -ts_adapt_type none -pc_type lu -ts_dt {{0.005 0.001 0.0005}separate output} -ts_max_steps {{100 500 1000}separate output} -ts_eimex_row_col {{1,1 2,1 3,1 2,2 3,2 3,3}separate output} 235c4762a1bSJed Brown requires: x 236c4762a1bSJed Brown 237c4762a1bSJed Brown TEST*/ 238