1c4762a1bSJed Brown static char help[] = "Solves the motion of spring.\n\
2c4762a1bSJed Brown Input parameters include:\n";
3c4762a1bSJed Brown
4c4762a1bSJed Brown /* ------------------------------------------------------------------------
5c4762a1bSJed Brown
6c4762a1bSJed Brown This program solves the motion of spring by Hooke's law
7c4762a1bSJed Brown x' = f(t,v) = v
8c4762a1bSJed Brown v' = g(t,x) = -omega^2*x
9c4762a1bSJed Brown on the interval 0 <= t <= 0.1, with the initial conditions
10c4762a1bSJed Brown x(0) = 0.2, x'(0) = v(0) = 0,
11c4762a1bSJed Brown and
12c4762a1bSJed Brown omega = 64.
13c4762a1bSJed Brown The exact solution is
14c4762a1bSJed Brown x(t) = A*sin(t*omega) + B*cos(t*omega)
15c4762a1bSJed Brown where A and B are constants that can be determined from the initial conditions.
16c4762a1bSJed Brown In this case, B=0.2, A=0.
17c4762a1bSJed Brown
18c4762a1bSJed Brown Notes:
19c4762a1bSJed Brown This code demonstrates the TS solver interface to solve a separable Hamiltonian
20c4762a1bSJed Brown system, which can be split into two subsystems involving two coupling components,
21c4762a1bSJed Brown named generailized momentum and generailized position respectively.
22c4762a1bSJed Brown Using a symplectic intergrator can preserve energy
23c4762a1bSJed Brown E = (v^2+omega^2*x^2-omega^2*h*v*x)/2
24c4762a1bSJed Brown ------------------------------------------------------------------------- */
25c4762a1bSJed Brown
26c4762a1bSJed Brown #include <petscts.h>
27c4762a1bSJed Brown #include <petscvec.h>
28c4762a1bSJed Brown
29c4762a1bSJed Brown typedef struct _n_User *User;
30c4762a1bSJed Brown struct _n_User {
31c4762a1bSJed Brown PetscReal omega;
32c4762a1bSJed Brown PetscInt nts; /* print the energy at each nts time steps */
33c4762a1bSJed Brown };
34c4762a1bSJed Brown
35c4762a1bSJed Brown /*
36c4762a1bSJed Brown User-defined routines.
37c4762a1bSJed Brown The first RHS function provides f(t,x), the residual for the generalized momentum,
38c4762a1bSJed Brown and the second one provides g(t,v), the residual for the generalized position.
39c4762a1bSJed Brown */
RHSFunction2(TS ts,PetscReal t,Vec X,Vec Vres,PetscCtx ctx)40*2a8381b2SBarry Smith static PetscErrorCode RHSFunction2(TS ts, PetscReal t, Vec X, Vec Vres, PetscCtx ctx)
41d71ae5a4SJacob Faibussowitsch {
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));
523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
53c4762a1bSJed Brown }
54c4762a1bSJed Brown
RHSFunction1(TS ts,PetscReal t,Vec V,Vec Xres,PetscCtx ctx)55*2a8381b2SBarry Smith static PetscErrorCode RHSFunction1(TS ts, PetscReal t, Vec V, Vec Xres, PetscCtx ctx)
56d71ae5a4SJacob Faibussowitsch {
57c4762a1bSJed Brown const PetscScalar *v;
58c4762a1bSJed Brown PetscScalar *xres;
59c4762a1bSJed Brown
60c4762a1bSJed Brown PetscFunctionBeginUser;
619566063dSJacob Faibussowitsch PetscCall(VecGetArray(Xres, &xres));
629566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(V, &v));
63c4762a1bSJed Brown xres[0] = v[0];
649566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(V, &v));
659566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Xres, &xres));
663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
67c4762a1bSJed Brown }
68c4762a1bSJed Brown
RHSFunction(TS ts,PetscReal t,Vec U,Vec R,PetscCtx ctx)69*2a8381b2SBarry Smith static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec U, Vec R, PetscCtx ctx)
70d71ae5a4SJacob Faibussowitsch {
71c4762a1bSJed Brown User user = (User)ctx;
72c4762a1bSJed Brown const PetscScalar *u;
73c4762a1bSJed Brown PetscScalar *r;
74c4762a1bSJed Brown
75c4762a1bSJed Brown PetscFunctionBeginUser;
769566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u));
779566063dSJacob Faibussowitsch PetscCall(VecGetArray(R, &r));
78c4762a1bSJed Brown r[0] = u[1];
79c4762a1bSJed Brown r[1] = -user->omega * user->omega * u[0];
809566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u));
819566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(R, &r));
823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
83c4762a1bSJed Brown }
84c4762a1bSJed Brown
85c4762a1bSJed Brown /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
Monitor(TS ts,PetscInt step,PetscReal t,Vec U,PetscCtx ctx)86*2a8381b2SBarry Smith static PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal t, Vec U, PetscCtx ctx)
87d71ae5a4SJacob Faibussowitsch {
88c4762a1bSJed Brown const PetscScalar *u;
89c4762a1bSJed Brown PetscReal dt;
90c4762a1bSJed Brown PetscScalar energy, menergy;
91c4762a1bSJed Brown User user = (User)ctx;
92c4762a1bSJed Brown
93c4762a1bSJed Brown PetscFunctionBeginUser;
94c4762a1bSJed Brown if (step % user->nts == 0) {
959566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts, &dt));
969566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u));
97c4762a1bSJed Brown menergy = (u[1] * u[1] + user->omega * user->omega * u[0] * u[0] - user->omega * user->omega * dt * u[0] * u[1]) / 2.;
98c4762a1bSJed Brown energy = (u[1] * u[1] + user->omega * user->omega * u[0] * u[0]) / 2.;
9963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "At time %.6lf, Energy = %8g, Modified Energy = %8g\n", (double)t, (double)energy, (double)menergy));
1009566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u));
101c4762a1bSJed Brown }
1023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
103c4762a1bSJed Brown }
104c4762a1bSJed Brown
main(int argc,char ** argv)105d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
106d71ae5a4SJacob Faibussowitsch {
107c4762a1bSJed Brown TS ts; /* nonlinear solver */
108c4762a1bSJed Brown Vec U; /* solution, residual vectors */
109c4762a1bSJed Brown IS is1, is2;
110c4762a1bSJed Brown PetscInt nindices[1];
111c4762a1bSJed Brown PetscReal ftime = 0.1;
112c4762a1bSJed Brown PetscBool monitor = PETSC_FALSE;
113c4762a1bSJed Brown PetscScalar *u_ptr;
114c4762a1bSJed Brown PetscMPIInt size;
115c4762a1bSJed Brown struct _n_User user;
116c4762a1bSJed Brown
117c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
118c4762a1bSJed Brown Initialize program
119c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
120327415f7SBarry Smith PetscFunctionBeginUser;
1219566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help));
1229566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
1233c633725SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
124c4762a1bSJed Brown
125c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
126c4762a1bSJed Brown Set runtime options
127c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
128c4762a1bSJed Brown user.omega = 64.;
129c4762a1bSJed Brown user.nts = 100;
1309566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-monitor", &monitor, NULL));
131d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Physical parameters", NULL);
132f3fa974cSJacob Faibussowitsch PetscCall(PetscOptionsReal("-omega", "parameter", "<64>", user.omega, &user.omega, NULL));
133f3fa974cSJacob Faibussowitsch PetscCall(PetscOptionsInt("-next_output", "time steps for next output point", "<100>", user.nts, &user.nts, NULL));
134d0609cedSBarry Smith PetscOptionsEnd();
135c4762a1bSJed Brown
136c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
137c4762a1bSJed Brown Create necessary matrix and vectors, solve same ODE on every process
138c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1399566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, 2, &U));
140c4762a1bSJed Brown nindices[0] = 0;
1419566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 1, nindices, PETSC_COPY_VALUES, &is1));
142c4762a1bSJed Brown nindices[0] = 1;
1439566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 1, nindices, PETSC_COPY_VALUES, &is2));
144c4762a1bSJed Brown
145c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
146c4762a1bSJed Brown Create timestepping solver context
147c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1489566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
1499566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSBASICSYMPLECTIC));
1509566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetIS(ts, "position", is1));
1519566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetIS(ts, "momentum", is2));
1529566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetRHSFunction(ts, "position", NULL, RHSFunction1, &user));
1539566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetRHSFunction(ts, "momentum", NULL, RHSFunction2, &user));
1549566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, &user));
155c4762a1bSJed Brown
1569566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, ftime));
1579566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, 0.0001));
1589566063dSJacob Faibussowitsch PetscCall(TSSetMaxSteps(ts, 1000));
1599566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
16048a46eb9SPierre Jolivet if (monitor) PetscCall(TSMonitorSet(ts, Monitor, &user, NULL));
161c4762a1bSJed Brown
162c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
163c4762a1bSJed Brown Set initial conditions
164c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1659566063dSJacob Faibussowitsch PetscCall(VecGetArray(U, &u_ptr));
166c4762a1bSJed Brown u_ptr[0] = 0.2;
167c4762a1bSJed Brown u_ptr[1] = 0.0;
1689566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U, &u_ptr));
169c4762a1bSJed Brown
1709566063dSJacob Faibussowitsch PetscCall(TSSetTime(ts, 0.0));
171c4762a1bSJed Brown
172c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
173c4762a1bSJed Brown Set runtime options
174c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1759566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts));
176c4762a1bSJed Brown
177c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
178c4762a1bSJed Brown Solve nonlinear system
179c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1809566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, U));
1819566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts, &ftime));
1829566063dSJacob Faibussowitsch PetscCall(VecView(U, PETSC_VIEWER_STDOUT_WORLD));
183c4762a1bSJed Brown
18463a3b9bcSJacob 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))));
185c4762a1bSJed Brown
186c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
187c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they
188c4762a1bSJed Brown are no longer needed.
189c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1909566063dSJacob Faibussowitsch PetscCall(VecDestroy(&U));
1919566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts));
1929566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is1));
1939566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is2));
1949566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
195b122ec5aSJacob Faibussowitsch return 0;
196c4762a1bSJed Brown }
197c4762a1bSJed Brown
198c4762a1bSJed Brown /*TEST
199c4762a1bSJed Brown build:
200c4762a1bSJed Brown requires: !single !complex
201c4762a1bSJed Brown
202c4762a1bSJed Brown test:
203c4762a1bSJed Brown args: -ts_basicsymplectic_type 1 -monitor
204c4762a1bSJed Brown
205c4762a1bSJed Brown test:
206c4762a1bSJed Brown suffix: 2
207c4762a1bSJed Brown args: -ts_basicsymplectic_type 2 -monitor
208c4762a1bSJed Brown
209c4762a1bSJed Brown TEST*/
210