xref: /petsc/src/ts/tutorials/ex49.c (revision 3ba1676111f5c958fe6c2729b46ca4d523958bb3)
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) = -6.6e-01,
12c4762a1bSJed Brown    and
13c4762a1bSJed Brown        mu = 10^6.
14c4762a1bSJed Brown    This is a nonlinear equation.
15c4762a1bSJed Brown 
16c4762a1bSJed Brown    This is a copy and modification of ex20.c to exactly match a test
17c4762a1bSJed Brown    problem that comes with the Radau5 integrator package.
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 
29d71ae5a4SJacob Faibussowitsch static PetscErrorCode IFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, void *ctx)
30d71ae5a4SJacob Faibussowitsch {
31c4762a1bSJed Brown   User               user = (User)ctx;
32c4762a1bSJed Brown   const PetscScalar *x, *xdot;
33c4762a1bSJed Brown   PetscScalar       *f;
34c4762a1bSJed Brown 
35c4762a1bSJed Brown   PetscFunctionBeginUser;
369566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
379566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xdot, &xdot));
389566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
39c4762a1bSJed Brown   f[0] = xdot[0] - x[1];
40c4762a1bSJed Brown   f[1] = xdot[1] - user->mu * ((1.0 - x[0] * x[0]) * x[1] - x[0]);
419566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
429566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xdot, &xdot));
439566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
44*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
45c4762a1bSJed Brown }
46c4762a1bSJed Brown 
47d71ae5a4SJacob Faibussowitsch static PetscErrorCode IJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat A, Mat B, void *ctx)
48d71ae5a4SJacob Faibussowitsch {
49c4762a1bSJed Brown   User               user     = (User)ctx;
50c4762a1bSJed Brown   PetscInt           rowcol[] = {0, 1};
51c4762a1bSJed Brown   const PetscScalar *x;
52c4762a1bSJed Brown   PetscScalar        J[2][2];
53c4762a1bSJed Brown 
54c4762a1bSJed Brown   PetscFunctionBeginUser;
559566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
569371c9d4SSatish Balay   J[0][0] = a;
579371c9d4SSatish Balay   J[0][1] = -1.0;
589371c9d4SSatish Balay   J[1][0] = user->mu * (1.0 + 2.0 * x[0] * x[1]);
599371c9d4SSatish Balay   J[1][1] = a - user->mu * (1.0 - x[0] * x[0]);
609566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES));
619566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
62c4762a1bSJed Brown 
639566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
649566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
65c4762a1bSJed Brown   if (A != B) {
669566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
679566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
68c4762a1bSJed Brown   }
69*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
70c4762a1bSJed Brown }
71c4762a1bSJed Brown 
72d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
73d71ae5a4SJacob Faibussowitsch {
74c4762a1bSJed Brown   TS             ts; /* nonlinear solver */
75c4762a1bSJed Brown   Vec            x;  /* solution, residual vectors */
76c4762a1bSJed Brown   Mat            A;  /* Jacobian matrix */
77c4762a1bSJed Brown   PetscInt       steps;
78c4762a1bSJed Brown   PetscReal      ftime = 2;
79c4762a1bSJed Brown   PetscScalar   *x_ptr;
80c4762a1bSJed Brown   PetscMPIInt    size;
81c4762a1bSJed Brown   struct _n_User user;
82c4762a1bSJed Brown 
83c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
84c4762a1bSJed Brown      Initialize program
85c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
86327415f7SBarry Smith   PetscFunctionBeginUser;
879566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
889566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
893c633725SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
90c4762a1bSJed Brown 
91c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
92c4762a1bSJed Brown     Set runtime options
93c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
94c4762a1bSJed Brown   user.next_output = 0.0;
95c4762a1bSJed Brown   user.mu          = 1.0e6;
96d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Physical parameters", NULL);
979566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mu", "Stiffness parameter", "<1.0e6>", user.mu, &user.mu, NULL));
98d0609cedSBarry Smith   PetscOptionsEnd();
99c4762a1bSJed Brown 
100c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
101c4762a1bSJed Brown     Create necessary matrix and vectors, solve same ODE on every process
102c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1039566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
1049566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, 2, 2));
1059566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
1069566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
107c4762a1bSJed Brown 
1089566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &x, NULL));
109c4762a1bSJed Brown 
110c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
111c4762a1bSJed Brown      Create timestepping solver context
112c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1139566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
1149566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSBEULER));
1159566063dSJacob Faibussowitsch   PetscCall(TSSetIFunction(ts, NULL, IFunction, &user));
1169566063dSJacob Faibussowitsch   PetscCall(TSSetIJacobian(ts, A, A, IJacobian, &user));
117c4762a1bSJed Brown 
1189566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, ftime));
1199566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
1209566063dSJacob Faibussowitsch   PetscCall(TSSetTolerances(ts, 1.e-4, NULL, 1.e-4, NULL));
121c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
122c4762a1bSJed Brown      Set initial conditions
123c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1249566063dSJacob Faibussowitsch   PetscCall(VecGetArray(x, &x_ptr));
1259371c9d4SSatish Balay   x_ptr[0] = 2.0;
1269371c9d4SSatish Balay   x_ptr[1] = -6.6e-01;
1279566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(x, &x_ptr));
1289566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, .000001));
129c4762a1bSJed Brown 
130c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
131c4762a1bSJed Brown      Set runtime options
132c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1339566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
134c4762a1bSJed Brown 
135c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
136c4762a1bSJed Brown      Solve nonlinear system
137c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1389566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, x));
1399566063dSJacob Faibussowitsch   PetscCall(TSGetSolveTime(ts, &ftime));
1409566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts, &steps));
14163a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "steps %" PetscInt_FMT ", ftime %g\n", steps, (double)ftime));
1429566063dSJacob Faibussowitsch   PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));
143c4762a1bSJed Brown 
144c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
145c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
146c4762a1bSJed Brown      are no longer needed.
147c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1489566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
1499566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
1509566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
151c4762a1bSJed Brown 
1529566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
153d0609cedSBarry Smith   return (0);
154c4762a1bSJed Brown }
155c4762a1bSJed Brown 
156c4762a1bSJed Brown /*TEST
157c4762a1bSJed Brown 
158c4762a1bSJed Brown     build:
159dfd57a17SPierre Jolivet       requires: double !complex !defined(PETSC_USE_64BIT_INDICES) radau5
160c4762a1bSJed Brown 
161c4762a1bSJed Brown     test:
162c4762a1bSJed Brown       args: -ts_monitor_solution -ts_type radau5
163c4762a1bSJed Brown 
164c4762a1bSJed Brown TEST*/
165