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 Concepts: TS^Separable Hamiltonian problems 7c4762a1bSJed Brown Concepts: TS^Symplectic intergrators 8c4762a1bSJed Brown Processors: 1 9c4762a1bSJed Brown */ 10c4762a1bSJed Brown /* ------------------------------------------------------------------------ 11c4762a1bSJed Brown 12c4762a1bSJed Brown This program solves the motion of spring by Hooke's law 13c4762a1bSJed Brown x' = f(t,v) = v 14c4762a1bSJed Brown v' = g(t,x) = -omega^2*x 15c4762a1bSJed Brown on the interval 0 <= t <= 0.1, with the initial conditions 16c4762a1bSJed Brown x(0) = 0.2, x'(0) = v(0) = 0, 17c4762a1bSJed Brown and 18c4762a1bSJed Brown omega = 64. 19c4762a1bSJed Brown The exact solution is 20c4762a1bSJed Brown x(t) = A*sin(t*omega) + B*cos(t*omega) 21c4762a1bSJed Brown where A and B are constants that can be determined from the initial conditions. 22c4762a1bSJed Brown In this case, B=0.2, A=0. 23c4762a1bSJed Brown 24c4762a1bSJed Brown Notes: 25c4762a1bSJed Brown This code demonstrates the TS solver interface to solve a separable Hamiltonian 26c4762a1bSJed Brown system, which can be split into two subsystems involving two coupling components, 27c4762a1bSJed Brown named generailized momentum and generailized position respectively. 28c4762a1bSJed Brown Using a symplectic intergrator can preserve energy 29c4762a1bSJed Brown E = (v^2+omega^2*x^2-omega^2*h*v*x)/2 30c4762a1bSJed Brown ------------------------------------------------------------------------- */ 31c4762a1bSJed Brown 32c4762a1bSJed Brown #include <petscts.h> 33c4762a1bSJed Brown #include <petscvec.h> 34c4762a1bSJed Brown 35c4762a1bSJed Brown typedef struct _n_User *User; 36c4762a1bSJed Brown struct _n_User { 37c4762a1bSJed Brown PetscReal omega; 38c4762a1bSJed Brown PetscInt nts; /* print the energy at each nts time steps */ 39c4762a1bSJed Brown }; 40c4762a1bSJed Brown 41c4762a1bSJed Brown /* 42c4762a1bSJed Brown User-defined routines. 43c4762a1bSJed Brown The first RHS function provides f(t,x), the residual for the generalized momentum, 44c4762a1bSJed Brown and the second one provides g(t,v), the residual for the generalized position. 45c4762a1bSJed Brown */ 46c4762a1bSJed Brown static PetscErrorCode RHSFunction2(TS ts,PetscReal t,Vec X,Vec Vres,void *ctx) 47c4762a1bSJed Brown { 48c4762a1bSJed Brown User user = (User)ctx; 49c4762a1bSJed Brown const PetscScalar *x; 50c4762a1bSJed Brown PetscScalar *vres; 51c4762a1bSJed Brown 52c4762a1bSJed Brown PetscFunctionBeginUser; 53*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&x)); 54*9566063dSJacob Faibussowitsch PetscCall(VecGetArray(Vres,&vres)); 55c4762a1bSJed Brown vres[0] = -user->omega*user->omega*x[0]; 56*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Vres,&vres)); 57*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&x)); 58c4762a1bSJed Brown PetscFunctionReturn(0); 59c4762a1bSJed Brown } 60c4762a1bSJed Brown 61c4762a1bSJed Brown static PetscErrorCode RHSFunction1(TS ts,PetscReal t,Vec V,Vec Xres,void *ctx) 62c4762a1bSJed Brown { 63c4762a1bSJed Brown const PetscScalar *v; 64c4762a1bSJed Brown PetscScalar *xres; 65c4762a1bSJed Brown 66c4762a1bSJed Brown PetscFunctionBeginUser; 67*9566063dSJacob Faibussowitsch PetscCall(VecGetArray(Xres,&xres)); 68*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(V,&v)); 69c4762a1bSJed Brown xres[0] = v[0]; 70*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(V,&v)); 71*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(Xres,&xres)); 72c4762a1bSJed Brown PetscFunctionReturn(0); 73c4762a1bSJed Brown } 74c4762a1bSJed Brown 75c4762a1bSJed Brown static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec U,Vec R,void *ctx) 76c4762a1bSJed Brown { 77c4762a1bSJed Brown User user = (User)ctx; 78c4762a1bSJed Brown const PetscScalar *u; 79c4762a1bSJed Brown PetscScalar *r; 80c4762a1bSJed Brown 81c4762a1bSJed Brown PetscFunctionBeginUser; 82*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U,&u)); 83*9566063dSJacob Faibussowitsch PetscCall(VecGetArray(R,&r)); 84c4762a1bSJed Brown r[0] = u[1]; 85c4762a1bSJed Brown r[1] = -user->omega*user->omega*u[0]; 86*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U,&u)); 87*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(R,&r)); 88c4762a1bSJed Brown PetscFunctionReturn(0); 89c4762a1bSJed Brown } 90c4762a1bSJed Brown 91c4762a1bSJed Brown /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */ 92c4762a1bSJed Brown static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec U,void *ctx) 93c4762a1bSJed Brown { 94c4762a1bSJed Brown const PetscScalar *u; 95c4762a1bSJed Brown PetscReal dt; 96c4762a1bSJed Brown PetscScalar energy,menergy; 97c4762a1bSJed Brown User user = (User)ctx; 98c4762a1bSJed Brown 99c4762a1bSJed Brown PetscFunctionBeginUser; 100c4762a1bSJed Brown if (step%user->nts == 0) { 101*9566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts,&dt)); 102*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U,&u)); 103c4762a1bSJed Brown menergy = (u[1]*u[1]+user->omega*user->omega*u[0]*u[0]-user->omega*user->omega*dt*u[0]*u[1])/2.; 104c4762a1bSJed Brown energy = (u[1]*u[1]+user->omega*user->omega*u[0]*u[0])/2.; 105*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"At time %.6lf, Energy = %8g, Modified Energy = %8g\n",t,(double)energy,(double)menergy)); 106*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U,&u)); 107c4762a1bSJed Brown } 108c4762a1bSJed Brown PetscFunctionReturn(0); 109c4762a1bSJed Brown } 110c4762a1bSJed Brown 111c4762a1bSJed Brown int main(int argc,char **argv) 112c4762a1bSJed Brown { 113c4762a1bSJed Brown TS ts; /* nonlinear solver */ 114c4762a1bSJed Brown Vec U; /* solution, residual vectors */ 115c4762a1bSJed Brown IS is1,is2; 116c4762a1bSJed Brown PetscInt nindices[1]; 117c4762a1bSJed Brown PetscReal ftime = 0.1; 118c4762a1bSJed Brown PetscBool monitor = PETSC_FALSE; 119c4762a1bSJed Brown PetscScalar *u_ptr; 120c4762a1bSJed Brown PetscMPIInt size; 121c4762a1bSJed Brown struct _n_User user; 122c4762a1bSJed Brown PetscErrorCode ierr; 123c4762a1bSJed Brown 124c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 125c4762a1bSJed Brown Initialize program 126c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 127*9566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,NULL,help)); 128*9566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 1293c633725SBarry Smith PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 130c4762a1bSJed Brown 131c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 132c4762a1bSJed Brown Set runtime options 133c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 134c4762a1bSJed Brown user.omega = 64.; 135c4762a1bSJed Brown user.nts = 100; 136*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL)); 137*9566063dSJacob Faibussowitsch ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Physical parameters",NULL);PetscCall(ierr); 138*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-omega","parameter","<64>",user.omega,&user.omega,PETSC_NULL)); 139*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-next_output","time steps for next output point","<100>",user.nts,&user.nts,PETSC_NULL)); 140*9566063dSJacob Faibussowitsch ierr = PetscOptionsEnd();PetscCall(ierr); 141c4762a1bSJed Brown 142c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 143c4762a1bSJed Brown Create necessary matrix and vectors, solve same ODE on every process 144c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 145*9566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,2,&U)); 146c4762a1bSJed Brown nindices[0] = 0; 147*9566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,1,nindices,PETSC_COPY_VALUES,&is1)); 148c4762a1bSJed Brown nindices[0] = 1; 149*9566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,1,nindices,PETSC_COPY_VALUES,&is2)); 150c4762a1bSJed Brown 151c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 152c4762a1bSJed Brown Create timestepping solver context 153c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 154*9566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD,&ts)); 155*9566063dSJacob Faibussowitsch PetscCall(TSSetType(ts,TSBASICSYMPLECTIC)); 156*9566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetIS(ts,"position",is1)); 157*9566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetIS(ts,"momentum",is2)); 158*9566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetRHSFunction(ts,"position",NULL,RHSFunction1,&user)); 159*9566063dSJacob Faibussowitsch PetscCall(TSRHSSplitSetRHSFunction(ts,"momentum",NULL,RHSFunction2,&user)); 160*9566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts,NULL,RHSFunction,&user)); 161c4762a1bSJed Brown 162*9566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts,ftime)); 163*9566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts,0.0001)); 164*9566063dSJacob Faibussowitsch PetscCall(TSSetMaxSteps(ts,1000)); 165*9566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP)); 166c4762a1bSJed Brown if (monitor) { 167*9566063dSJacob Faibussowitsch PetscCall(TSMonitorSet(ts,Monitor,&user,NULL)); 168c4762a1bSJed Brown } 169c4762a1bSJed Brown 170c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 171c4762a1bSJed Brown Set initial conditions 172c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 173*9566063dSJacob Faibussowitsch PetscCall(VecGetArray(U,&u_ptr)); 174c4762a1bSJed Brown u_ptr[0] = 0.2; 175c4762a1bSJed Brown u_ptr[1] = 0.0; 176*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U,&u_ptr)); 177c4762a1bSJed Brown 178*9566063dSJacob Faibussowitsch PetscCall(TSSetTime(ts,0.0)); 179c4762a1bSJed Brown 180c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 181c4762a1bSJed Brown Set runtime options 182c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 183*9566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 184c4762a1bSJed Brown 185c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 186c4762a1bSJed Brown Solve nonlinear system 187c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 188*9566063dSJacob Faibussowitsch PetscCall(TSSolve(ts,U)); 189*9566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts,&ftime)); 190*9566063dSJacob Faibussowitsch PetscCall(VecView(U,PETSC_VIEWER_STDOUT_WORLD)); 191c4762a1bSJed Brown 192*9566063dSJacob 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))); 193c4762a1bSJed Brown 194c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 195c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 196c4762a1bSJed Brown are no longer needed. 197c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 198*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&U)); 199*9566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 200*9566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is1)); 201*9566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is2)); 202*9566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 203b122ec5aSJacob Faibussowitsch return 0; 204c4762a1bSJed Brown } 205c4762a1bSJed Brown 206c4762a1bSJed Brown /*TEST 207c4762a1bSJed Brown build: 208c4762a1bSJed Brown requires: !single !complex 209c4762a1bSJed Brown 210c4762a1bSJed Brown test: 211c4762a1bSJed Brown args: -ts_basicsymplectic_type 1 -monitor 212c4762a1bSJed Brown 213c4762a1bSJed Brown test: 214c4762a1bSJed Brown suffix: 2 215c4762a1bSJed Brown args: -ts_basicsymplectic_type 2 -monitor 216c4762a1bSJed Brown 217c4762a1bSJed Brown TEST*/ 218