xref: /petsc/src/ts/tutorials/hamiltonian/ex1.c (revision d0609ced746bc51b019815ca91d747429db24893)
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 */
41c4762a1bSJed Brown static PetscErrorCode RHSFunction2(TS ts,PetscReal t,Vec X,Vec Vres,void *ctx)
42c4762a1bSJed Brown {
43c4762a1bSJed Brown   User              user = (User)ctx;
44c4762a1bSJed Brown   const PetscScalar *x;
45c4762a1bSJed Brown   PetscScalar       *vres;
46c4762a1bSJed Brown 
47c4762a1bSJed Brown   PetscFunctionBeginUser;
489566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X,&x));
499566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Vres,&vres));
50c4762a1bSJed Brown   vres[0] = -user->omega*user->omega*x[0];
519566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Vres,&vres));
529566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X,&x));
53c4762a1bSJed Brown   PetscFunctionReturn(0);
54c4762a1bSJed Brown }
55c4762a1bSJed Brown 
56c4762a1bSJed Brown static PetscErrorCode RHSFunction1(TS ts,PetscReal t,Vec V,Vec Xres,void *ctx)
57c4762a1bSJed Brown {
58c4762a1bSJed Brown   const PetscScalar *v;
59c4762a1bSJed Brown   PetscScalar       *xres;
60c4762a1bSJed Brown 
61c4762a1bSJed Brown   PetscFunctionBeginUser;
629566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Xres,&xres));
639566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(V,&v));
64c4762a1bSJed Brown   xres[0] = v[0];
659566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(V,&v));
669566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Xres,&xres));
67c4762a1bSJed Brown   PetscFunctionReturn(0);
68c4762a1bSJed Brown }
69c4762a1bSJed Brown 
70c4762a1bSJed Brown static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec U,Vec R,void *ctx)
71c4762a1bSJed Brown {
72c4762a1bSJed Brown   User              user = (User)ctx;
73c4762a1bSJed Brown   const PetscScalar *u;
74c4762a1bSJed Brown   PetscScalar       *r;
75c4762a1bSJed Brown 
76c4762a1bSJed Brown   PetscFunctionBeginUser;
779566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U,&u));
789566063dSJacob Faibussowitsch   PetscCall(VecGetArray(R,&r));
79c4762a1bSJed Brown   r[0] = u[1];
80c4762a1bSJed Brown   r[1] = -user->omega*user->omega*u[0];
819566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U,&u));
829566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(R,&r));
83c4762a1bSJed Brown   PetscFunctionReturn(0);
84c4762a1bSJed Brown }
85c4762a1bSJed Brown 
86c4762a1bSJed Brown /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
87c4762a1bSJed Brown static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec U,void *ctx)
88c4762a1bSJed Brown {
89c4762a1bSJed Brown   const PetscScalar *u;
90c4762a1bSJed Brown   PetscReal         dt;
91c4762a1bSJed Brown   PetscScalar       energy,menergy;
92c4762a1bSJed Brown   User              user = (User)ctx;
93c4762a1bSJed Brown 
94c4762a1bSJed Brown   PetscFunctionBeginUser;
95c4762a1bSJed Brown   if (step%user->nts == 0) {
969566063dSJacob Faibussowitsch     PetscCall(TSGetTimeStep(ts,&dt));
979566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(U,&u));
98c4762a1bSJed Brown     menergy = (u[1]*u[1]+user->omega*user->omega*u[0]*u[0]-user->omega*user->omega*dt*u[0]*u[1])/2.;
99c4762a1bSJed Brown     energy = (u[1]*u[1]+user->omega*user->omega*u[0]*u[0])/2.;
1009566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"At time %.6lf, Energy = %8g, Modified Energy = %8g\n",t,(double)energy,(double)menergy));
1019566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(U,&u));
102c4762a1bSJed Brown   }
103c4762a1bSJed Brown   PetscFunctionReturn(0);
104c4762a1bSJed Brown }
105c4762a1bSJed Brown 
106c4762a1bSJed Brown int main(int argc,char **argv)
107c4762a1bSJed Brown {
108c4762a1bSJed Brown   TS             ts;            /* nonlinear solver */
109c4762a1bSJed Brown   Vec            U;             /* solution, residual vectors */
110c4762a1bSJed Brown   IS             is1,is2;
111c4762a1bSJed Brown   PetscInt       nindices[1];
112c4762a1bSJed Brown   PetscReal      ftime   = 0.1;
113c4762a1bSJed Brown   PetscBool      monitor = PETSC_FALSE;
114c4762a1bSJed Brown   PetscScalar    *u_ptr;
115c4762a1bSJed Brown   PetscMPIInt    size;
116c4762a1bSJed Brown   struct _n_User user;
117c4762a1bSJed Brown 
118c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
119c4762a1bSJed Brown      Initialize program
120c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
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));
131*d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Physical parameters",NULL);
1329566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-omega","parameter","<64>",user.omega,&user.omega,PETSC_NULL));
1339566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-next_output","time steps for next output point","<100>",user.nts,&user.nts,PETSC_NULL));
134*d0609cedSBarry 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));
160c4762a1bSJed Brown   if (monitor) {
1619566063dSJacob Faibussowitsch     PetscCall(TSMonitorSet(ts,Monitor,&user,NULL));
162c4762a1bSJed Brown   }
163c4762a1bSJed Brown 
164c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
165c4762a1bSJed Brown      Set initial conditions
166c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1679566063dSJacob Faibussowitsch   PetscCall(VecGetArray(U,&u_ptr));
168c4762a1bSJed Brown   u_ptr[0] = 0.2;
169c4762a1bSJed Brown   u_ptr[1] = 0.0;
1709566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(U,&u_ptr));
171c4762a1bSJed Brown 
1729566063dSJacob Faibussowitsch   PetscCall(TSSetTime(ts,0.0));
173c4762a1bSJed Brown 
174c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
175c4762a1bSJed Brown      Set runtime options
176c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1779566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
178c4762a1bSJed Brown 
179c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
180c4762a1bSJed Brown      Solve nonlinear system
181c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1829566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts,U));
1839566063dSJacob Faibussowitsch   PetscCall(TSGetSolveTime(ts,&ftime));
1849566063dSJacob Faibussowitsch   PetscCall(VecView(U,PETSC_VIEWER_STDOUT_WORLD));
185c4762a1bSJed Brown 
1869566063dSJacob 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)));
187c4762a1bSJed Brown 
188c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
189c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
190c4762a1bSJed Brown      are no longer needed.
191c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1929566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&U));
1939566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
1949566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is1));
1959566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is2));
1969566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
197b122ec5aSJacob Faibussowitsch   return 0;
198c4762a1bSJed Brown }
199c4762a1bSJed Brown 
200c4762a1bSJed Brown /*TEST
201c4762a1bSJed Brown    build:
202c4762a1bSJed Brown      requires: !single !complex
203c4762a1bSJed Brown 
204c4762a1bSJed Brown    test:
205c4762a1bSJed Brown      args: -ts_basicsymplectic_type 1 -monitor
206c4762a1bSJed Brown 
207c4762a1bSJed Brown    test:
208c4762a1bSJed Brown      suffix: 2
209c4762a1bSJed Brown      args: -ts_basicsymplectic_type 2 -monitor
210c4762a1bSJed Brown 
211c4762a1bSJed Brown TEST*/
212