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) = -6.6e-01, 11c4762a1bSJed Brown and 12c4762a1bSJed Brown mu = 10^6. 13c4762a1bSJed Brown This is a nonlinear equation. 14c4762a1bSJed Brown 15c4762a1bSJed Brown This is a copy and modification of ex20.c to exactly match a test 16c4762a1bSJed Brown problem that comes with the Radau5 integrator package. 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 28*2a8381b2SBarry Smith static PetscErrorCode IFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, PetscCtx ctx) 29d71ae5a4SJacob Faibussowitsch { 30c4762a1bSJed Brown User user = (User)ctx; 31c4762a1bSJed Brown const PetscScalar *x, *xdot; 32c4762a1bSJed Brown PetscScalar *f; 33c4762a1bSJed Brown 34c4762a1bSJed Brown PetscFunctionBeginUser; 359566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 369566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Xdot, &xdot)); 379566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 38c4762a1bSJed Brown f[0] = xdot[0] - x[1]; 39c4762a1bSJed Brown f[1] = xdot[1] - user->mu * ((1.0 - x[0] * x[0]) * x[1] - x[0]); 409566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 419566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Xdot, &xdot)); 429566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 44c4762a1bSJed Brown } 45c4762a1bSJed Brown 46*2a8381b2SBarry Smith static PetscErrorCode IJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat A, Mat B, PetscCtx ctx) 47d71ae5a4SJacob Faibussowitsch { 48c4762a1bSJed Brown User user = (User)ctx; 49c4762a1bSJed Brown PetscInt rowcol[] = {0, 1}; 50c4762a1bSJed Brown const PetscScalar *x; 51c4762a1bSJed Brown PetscScalar J[2][2]; 52c4762a1bSJed Brown 53c4762a1bSJed Brown PetscFunctionBeginUser; 549566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 559371c9d4SSatish Balay J[0][0] = a; 569371c9d4SSatish Balay J[0][1] = -1.0; 579371c9d4SSatish Balay J[1][0] = user->mu * (1.0 + 2.0 * x[0] * x[1]); 589371c9d4SSatish Balay J[1][1] = a - user->mu * (1.0 - x[0] * x[0]); 599566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES)); 609566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 61c4762a1bSJed Brown 629566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 639566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 64c4762a1bSJed Brown if (A != B) { 659566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 669566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 67c4762a1bSJed Brown } 683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 69c4762a1bSJed Brown } 70c4762a1bSJed Brown 71d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 72d71ae5a4SJacob Faibussowitsch { 73c4762a1bSJed Brown TS ts; /* nonlinear solver */ 74c4762a1bSJed Brown Vec x; /* solution, residual vectors */ 75c4762a1bSJed Brown Mat A; /* Jacobian matrix */ 76c4762a1bSJed Brown PetscInt steps; 77c4762a1bSJed Brown PetscReal ftime = 2; 78c4762a1bSJed Brown PetscScalar *x_ptr; 79c4762a1bSJed Brown PetscMPIInt size; 80c4762a1bSJed Brown struct _n_User user; 81c4762a1bSJed Brown 82c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 83c4762a1bSJed Brown Initialize program 84c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 85327415f7SBarry Smith PetscFunctionBeginUser; 869566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 879566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 883c633725SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!"); 89c4762a1bSJed Brown 90c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 91c4762a1bSJed Brown Set runtime options 92c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 93c4762a1bSJed Brown user.next_output = 0.0; 94c4762a1bSJed Brown user.mu = 1.0e6; 95d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Physical parameters", NULL); 969566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-mu", "Stiffness parameter", "<1.0e6>", user.mu, &user.mu, NULL)); 97d0609cedSBarry Smith PetscOptionsEnd(); 98c4762a1bSJed Brown 99c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 100c4762a1bSJed Brown Create necessary matrix and vectors, solve same ODE on every process 101c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1029566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 1039566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, 2, 2)); 1049566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 1059566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 106c4762a1bSJed Brown 1079566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &x, NULL)); 108c4762a1bSJed Brown 109c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 110c4762a1bSJed Brown Create timestepping solver context 111c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1129566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 1139566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSBEULER)); 1149566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts, NULL, IFunction, &user)); 1159566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts, A, A, IJacobian, &user)); 116c4762a1bSJed Brown 1179566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, ftime)); 1189566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER)); 1199566063dSJacob Faibussowitsch PetscCall(TSSetTolerances(ts, 1.e-4, NULL, 1.e-4, NULL)); 120c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 121c4762a1bSJed Brown Set initial conditions 122c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1239566063dSJacob Faibussowitsch PetscCall(VecGetArray(x, &x_ptr)); 1249371c9d4SSatish Balay x_ptr[0] = 2.0; 1259371c9d4SSatish Balay x_ptr[1] = -6.6e-01; 1269566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x, &x_ptr)); 1279566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, .000001)); 128c4762a1bSJed Brown 129c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 130c4762a1bSJed Brown Set runtime options 131c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1329566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 133c4762a1bSJed Brown 134c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 135c4762a1bSJed Brown Solve nonlinear system 136c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1379566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, x)); 1389566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts, &ftime)); 1399566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &steps)); 14063a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "steps %" PetscInt_FMT ", ftime %g\n", steps, (double)ftime)); 1419566063dSJacob Faibussowitsch PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD)); 142c4762a1bSJed Brown 143c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 144c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 145c4762a1bSJed Brown are no longer needed. 146c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1479566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1489566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1499566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 150c4762a1bSJed Brown 1519566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 1524ad8454bSPierre Jolivet return 0; 153c4762a1bSJed Brown } 154c4762a1bSJed Brown 155c4762a1bSJed Brown /*TEST 156c4762a1bSJed Brown 157c4762a1bSJed Brown build: 158dfd57a17SPierre Jolivet requires: double !complex !defined(PETSC_USE_64BIT_INDICES) radau5 159c4762a1bSJed Brown 160c4762a1bSJed Brown test: 161c4762a1bSJed Brown args: -ts_monitor_solution -ts_type radau5 162c4762a1bSJed Brown 163c4762a1bSJed Brown TEST*/ 164