xref: /petsc/src/ts/tutorials/hamiltonian/ex1.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Solves the motion of spring.\n\
3c4762a1bSJed Brown Input parameters include:\n";
4c4762a1bSJed Brown 
5c4762a1bSJed Brown /* ------------------------------------------------------------------------
6c4762a1bSJed Brown 
7c4762a1bSJed Brown   This program solves the motion of spring by Hooke's law
8c4762a1bSJed Brown   x' = f(t,v) = v
9c4762a1bSJed Brown   v' = g(t,x) = -omega^2*x
10c4762a1bSJed Brown   on the interval 0 <= t <= 0.1, with the initial conditions
11c4762a1bSJed Brown     x(0) = 0.2, x'(0) = v(0) = 0,
12c4762a1bSJed Brown   and
13c4762a1bSJed Brown     omega = 64.
14c4762a1bSJed Brown   The exact solution is
15c4762a1bSJed Brown     x(t) = A*sin(t*omega) + B*cos(t*omega)
16c4762a1bSJed Brown   where A and B are constants that can be determined from the initial conditions.
17c4762a1bSJed Brown   In this case, B=0.2, A=0.
18c4762a1bSJed Brown 
19c4762a1bSJed Brown   Notes:
20c4762a1bSJed Brown   This code demonstrates the TS solver interface to solve a separable Hamiltonian
21c4762a1bSJed Brown   system, which can be split into two subsystems involving two coupling components,
22c4762a1bSJed Brown   named generailized momentum and generailized position respectively.
23c4762a1bSJed Brown   Using a symplectic intergrator can preserve energy
24c4762a1bSJed Brown   E = (v^2+omega^2*x^2-omega^2*h*v*x)/2
25c4762a1bSJed Brown   ------------------------------------------------------------------------- */
26c4762a1bSJed Brown 
27c4762a1bSJed Brown #include <petscts.h>
28c4762a1bSJed Brown #include <petscvec.h>
29c4762a1bSJed Brown 
30c4762a1bSJed Brown typedef struct _n_User *User;
31c4762a1bSJed Brown struct _n_User {
32c4762a1bSJed Brown   PetscReal omega;
33c4762a1bSJed Brown   PetscInt  nts; /* print the energy at each nts time steps */
34c4762a1bSJed Brown };
35c4762a1bSJed Brown 
36c4762a1bSJed Brown /*
37c4762a1bSJed Brown   User-defined routines.
38c4762a1bSJed Brown   The first RHS function provides f(t,x), the residual for the generalized momentum,
39c4762a1bSJed Brown   and the second one provides g(t,v), the residual for the generalized position.
40c4762a1bSJed Brown */
419371c9d4SSatish Balay static PetscErrorCode RHSFunction2(TS ts, PetscReal t, Vec X, Vec Vres, void *ctx) {
42c4762a1bSJed Brown   User               user = (User)ctx;
43c4762a1bSJed Brown   const PetscScalar *x;
44c4762a1bSJed Brown   PetscScalar       *vres;
45c4762a1bSJed Brown 
46c4762a1bSJed Brown   PetscFunctionBeginUser;
479566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
489566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Vres, &vres));
49c4762a1bSJed Brown   vres[0] = -user->omega * user->omega * x[0];
509566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Vres, &vres));
519566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
52c4762a1bSJed Brown   PetscFunctionReturn(0);
53c4762a1bSJed Brown }
54c4762a1bSJed Brown 
559371c9d4SSatish Balay static PetscErrorCode RHSFunction1(TS ts, PetscReal t, Vec V, Vec Xres, void *ctx) {
56c4762a1bSJed Brown   const PetscScalar *v;
57c4762a1bSJed Brown   PetscScalar       *xres;
58c4762a1bSJed Brown 
59c4762a1bSJed Brown   PetscFunctionBeginUser;
609566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Xres, &xres));
619566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(V, &v));
62c4762a1bSJed Brown   xres[0] = v[0];
639566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(V, &v));
649566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Xres, &xres));
65c4762a1bSJed Brown   PetscFunctionReturn(0);
66c4762a1bSJed Brown }
67c4762a1bSJed Brown 
689371c9d4SSatish Balay static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec U, Vec R, void *ctx) {
69c4762a1bSJed Brown   User               user = (User)ctx;
70c4762a1bSJed Brown   const PetscScalar *u;
71c4762a1bSJed Brown   PetscScalar       *r;
72c4762a1bSJed Brown 
73c4762a1bSJed Brown   PetscFunctionBeginUser;
749566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
759566063dSJacob Faibussowitsch   PetscCall(VecGetArray(R, &r));
76c4762a1bSJed Brown   r[0] = u[1];
77c4762a1bSJed Brown   r[1] = -user->omega * user->omega * u[0];
789566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
799566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(R, &r));
80c4762a1bSJed Brown   PetscFunctionReturn(0);
81c4762a1bSJed Brown }
82c4762a1bSJed Brown 
83c4762a1bSJed Brown /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
849371c9d4SSatish Balay static PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal t, Vec U, void *ctx) {
85c4762a1bSJed Brown   const PetscScalar *u;
86c4762a1bSJed Brown   PetscReal          dt;
87c4762a1bSJed Brown   PetscScalar        energy, menergy;
88c4762a1bSJed Brown   User               user = (User)ctx;
89c4762a1bSJed Brown 
90c4762a1bSJed Brown   PetscFunctionBeginUser;
91c4762a1bSJed Brown   if (step % user->nts == 0) {
929566063dSJacob Faibussowitsch     PetscCall(TSGetTimeStep(ts, &dt));
939566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(U, &u));
94c4762a1bSJed Brown     menergy = (u[1] * u[1] + user->omega * user->omega * u[0] * u[0] - user->omega * user->omega * dt * u[0] * u[1]) / 2.;
95c4762a1bSJed Brown     energy  = (u[1] * u[1] + user->omega * user->omega * u[0] * u[0]) / 2.;
9663a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "At time %.6lf, Energy = %8g, Modified Energy = %8g\n", (double)t, (double)energy, (double)menergy));
979566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(U, &u));
98c4762a1bSJed Brown   }
99c4762a1bSJed Brown   PetscFunctionReturn(0);
100c4762a1bSJed Brown }
101c4762a1bSJed Brown 
1029371c9d4SSatish Balay int main(int argc, char **argv) {
103c4762a1bSJed Brown   TS             ts; /* nonlinear solver */
104c4762a1bSJed Brown   Vec            U;  /* solution, residual vectors */
105c4762a1bSJed Brown   IS             is1, is2;
106c4762a1bSJed Brown   PetscInt       nindices[1];
107c4762a1bSJed Brown   PetscReal      ftime   = 0.1;
108c4762a1bSJed Brown   PetscBool      monitor = PETSC_FALSE;
109c4762a1bSJed Brown   PetscScalar   *u_ptr;
110c4762a1bSJed Brown   PetscMPIInt    size;
111c4762a1bSJed Brown   struct _n_User user;
112c4762a1bSJed Brown 
113c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
114c4762a1bSJed Brown      Initialize program
115c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
116327415f7SBarry Smith   PetscFunctionBeginUser;
1179566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
1189566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
1193c633725SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
120c4762a1bSJed Brown 
121c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
122c4762a1bSJed Brown     Set runtime options
123c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
124c4762a1bSJed Brown   user.omega = 64.;
125c4762a1bSJed Brown   user.nts   = 100;
1269566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-monitor", &monitor, NULL));
127d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Physical parameters", NULL);
1289566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-omega", "parameter", "<64>", user.omega, &user.omega, PETSC_NULL));
1299566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-next_output", "time steps for next output point", "<100>", user.nts, &user.nts, PETSC_NULL));
130d0609cedSBarry Smith   PetscOptionsEnd();
131c4762a1bSJed Brown 
132c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
133c4762a1bSJed Brown     Create necessary matrix and vectors, solve same ODE on every process
134c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1359566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_SELF, 2, &U));
136c4762a1bSJed Brown   nindices[0] = 0;
1379566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 1, nindices, PETSC_COPY_VALUES, &is1));
138c4762a1bSJed Brown   nindices[0] = 1;
1399566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 1, nindices, PETSC_COPY_VALUES, &is2));
140c4762a1bSJed Brown 
141c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
142c4762a1bSJed Brown      Create timestepping solver context
143c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1449566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
1459566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSBASICSYMPLECTIC));
1469566063dSJacob Faibussowitsch   PetscCall(TSRHSSplitSetIS(ts, "position", is1));
1479566063dSJacob Faibussowitsch   PetscCall(TSRHSSplitSetIS(ts, "momentum", is2));
1489566063dSJacob Faibussowitsch   PetscCall(TSRHSSplitSetRHSFunction(ts, "position", NULL, RHSFunction1, &user));
1499566063dSJacob Faibussowitsch   PetscCall(TSRHSSplitSetRHSFunction(ts, "momentum", NULL, RHSFunction2, &user));
1509566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, &user));
151c4762a1bSJed Brown 
1529566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, ftime));
1539566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, 0.0001));
1549566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSteps(ts, 1000));
1559566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
156*48a46eb9SPierre Jolivet   if (monitor) PetscCall(TSMonitorSet(ts, Monitor, &user, NULL));
157c4762a1bSJed Brown 
158c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
159c4762a1bSJed Brown      Set initial conditions
160c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1619566063dSJacob Faibussowitsch   PetscCall(VecGetArray(U, &u_ptr));
162c4762a1bSJed Brown   u_ptr[0] = 0.2;
163c4762a1bSJed Brown   u_ptr[1] = 0.0;
1649566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(U, &u_ptr));
165c4762a1bSJed Brown 
1669566063dSJacob Faibussowitsch   PetscCall(TSSetTime(ts, 0.0));
167c4762a1bSJed Brown 
168c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
169c4762a1bSJed Brown      Set runtime options
170c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1719566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
172c4762a1bSJed Brown 
173c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
174c4762a1bSJed Brown      Solve nonlinear system
175c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1769566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, U));
1779566063dSJacob Faibussowitsch   PetscCall(TSGetSolveTime(ts, &ftime));
1789566063dSJacob Faibussowitsch   PetscCall(VecView(U, PETSC_VIEWER_STDOUT_WORLD));
179c4762a1bSJed Brown 
18063a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "The exact solution at time %.6lf is [%g %g]\n", (double)ftime, (double)(0.2 * PetscCosReal(user.omega * ftime)), (double)(-0.2 * user.omega * PetscSinReal(user.omega * ftime))));
181c4762a1bSJed Brown 
182c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
183c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
184c4762a1bSJed Brown      are no longer needed.
185c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1869566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&U));
1879566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
1889566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is1));
1899566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is2));
1909566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
191b122ec5aSJacob Faibussowitsch   return 0;
192c4762a1bSJed Brown }
193c4762a1bSJed Brown 
194c4762a1bSJed Brown /*TEST
195c4762a1bSJed Brown    build:
196c4762a1bSJed Brown      requires: !single !complex
197c4762a1bSJed Brown 
198c4762a1bSJed Brown    test:
199c4762a1bSJed Brown      args: -ts_basicsymplectic_type 1 -monitor
200c4762a1bSJed Brown 
201c4762a1bSJed Brown    test:
202c4762a1bSJed Brown      suffix: 2
203c4762a1bSJed Brown      args: -ts_basicsymplectic_type 2 -monitor
204c4762a1bSJed Brown 
205c4762a1bSJed Brown TEST*/
206