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 PetscErrorCode ierr; 52c4762a1bSJed Brown 53c4762a1bSJed Brown PetscFunctionBeginUser; 54c4762a1bSJed Brown ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 55c4762a1bSJed Brown ierr = VecGetArray(Vres,&vres);CHKERRQ(ierr); 56c4762a1bSJed Brown vres[0] = -user->omega*user->omega*x[0]; 57c4762a1bSJed Brown ierr = VecRestoreArray(Vres,&vres);CHKERRQ(ierr); 58c4762a1bSJed Brown ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 59c4762a1bSJed Brown PetscFunctionReturn(0); 60c4762a1bSJed Brown } 61c4762a1bSJed Brown 62c4762a1bSJed Brown static PetscErrorCode RHSFunction1(TS ts,PetscReal t,Vec V,Vec Xres,void *ctx) 63c4762a1bSJed Brown { 64c4762a1bSJed Brown const PetscScalar *v; 65c4762a1bSJed Brown PetscScalar *xres; 66c4762a1bSJed Brown PetscErrorCode ierr; 67c4762a1bSJed Brown 68c4762a1bSJed Brown PetscFunctionBeginUser; 69c4762a1bSJed Brown ierr = VecGetArray(Xres,&xres);CHKERRQ(ierr); 70c4762a1bSJed Brown ierr = VecGetArrayRead(V,&v);CHKERRQ(ierr); 71c4762a1bSJed Brown xres[0] = v[0]; 72c4762a1bSJed Brown ierr = VecRestoreArrayRead(V,&v);CHKERRQ(ierr); 73c4762a1bSJed Brown ierr = VecRestoreArray(Xres,&xres);CHKERRQ(ierr); 74c4762a1bSJed Brown PetscFunctionReturn(0); 75c4762a1bSJed Brown } 76c4762a1bSJed Brown 77c4762a1bSJed Brown static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec U,Vec R,void *ctx) 78c4762a1bSJed Brown { 79c4762a1bSJed Brown User user = (User)ctx; 80c4762a1bSJed Brown const PetscScalar *u; 81c4762a1bSJed Brown PetscScalar *r; 82c4762a1bSJed Brown PetscErrorCode ierr; 83c4762a1bSJed Brown 84c4762a1bSJed Brown PetscFunctionBeginUser; 85c4762a1bSJed Brown ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 86c4762a1bSJed Brown ierr = VecGetArray(R,&r);CHKERRQ(ierr); 87c4762a1bSJed Brown r[0] = u[1]; 88c4762a1bSJed Brown r[1] = -user->omega*user->omega*u[0]; 89c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 90c4762a1bSJed Brown ierr = VecRestoreArray(R,&r);CHKERRQ(ierr); 91c4762a1bSJed Brown PetscFunctionReturn(0); 92c4762a1bSJed Brown } 93c4762a1bSJed Brown 94c4762a1bSJed Brown /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */ 95c4762a1bSJed Brown static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec U,void *ctx) 96c4762a1bSJed Brown { 97c4762a1bSJed Brown PetscErrorCode ierr; 98c4762a1bSJed Brown const PetscScalar *u; 99c4762a1bSJed Brown PetscReal dt; 100c4762a1bSJed Brown PetscScalar energy,menergy; 101c4762a1bSJed Brown User user = (User)ctx; 102c4762a1bSJed Brown 103c4762a1bSJed Brown PetscFunctionBeginUser; 104c4762a1bSJed Brown if (step%user->nts == 0) { 105c4762a1bSJed Brown ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr); 106c4762a1bSJed Brown ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 107c4762a1bSJed Brown menergy = (u[1]*u[1]+user->omega*user->omega*u[0]*u[0]-user->omega*user->omega*dt*u[0]*u[1])/2.; 108c4762a1bSJed Brown energy = (u[1]*u[1]+user->omega*user->omega*u[0]*u[0])/2.; 109c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"At time %.6lf, Energy = %8g, Modified Energy = %8g\n",t,(double)energy,(double)menergy);CHKERRQ(ierr); 110c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 111c4762a1bSJed Brown } 112c4762a1bSJed Brown PetscFunctionReturn(0); 113c4762a1bSJed Brown } 114c4762a1bSJed Brown 115c4762a1bSJed Brown int main(int argc,char **argv) 116c4762a1bSJed Brown { 117c4762a1bSJed Brown TS ts; /* nonlinear solver */ 118c4762a1bSJed Brown Vec U; /* solution, residual vectors */ 119c4762a1bSJed Brown IS is1,is2; 120c4762a1bSJed Brown PetscInt nindices[1]; 121c4762a1bSJed Brown PetscReal ftime = 0.1; 122c4762a1bSJed Brown PetscBool monitor = PETSC_FALSE; 123c4762a1bSJed Brown PetscScalar *u_ptr; 124c4762a1bSJed Brown PetscMPIInt size; 125c4762a1bSJed Brown struct _n_User user; 126c4762a1bSJed Brown PetscErrorCode ierr; 127c4762a1bSJed Brown 128c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 129c4762a1bSJed Brown Initialize program 130c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 131c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr; 132ffc4695bSBarry Smith ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 133*3c633725SBarry Smith PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 134c4762a1bSJed Brown 135c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 136c4762a1bSJed Brown Set runtime options 137c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 138c4762a1bSJed Brown user.omega = 64.; 139c4762a1bSJed Brown user.nts = 100; 140c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL);CHKERRQ(ierr); 141c4762a1bSJed Brown ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Physical parameters",NULL);CHKERRQ(ierr); 142c4762a1bSJed Brown ierr = PetscOptionsReal("-omega","parameter","<64>",user.omega,&user.omega,PETSC_NULL);CHKERRQ(ierr); 143c4762a1bSJed Brown ierr = PetscOptionsInt("-next_output","time steps for next output point","<100>",user.nts,&user.nts,PETSC_NULL);CHKERRQ(ierr); 144c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 145c4762a1bSJed Brown 146c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 147c4762a1bSJed Brown Create necessary matrix and vectors, solve same ODE on every process 148c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 149c4762a1bSJed Brown ierr = VecCreateSeq(PETSC_COMM_SELF,2,&U);CHKERRQ(ierr); 150c4762a1bSJed Brown nindices[0] = 0; 151c4762a1bSJed Brown ierr = ISCreateGeneral(PETSC_COMM_SELF,1,nindices,PETSC_COPY_VALUES,&is1);CHKERRQ(ierr); 152c4762a1bSJed Brown nindices[0] = 1; 153c4762a1bSJed Brown ierr = ISCreateGeneral(PETSC_COMM_SELF,1,nindices,PETSC_COPY_VALUES,&is2);CHKERRQ(ierr); 154c4762a1bSJed Brown 155c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 156c4762a1bSJed Brown Create timestepping solver context 157c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 158c4762a1bSJed Brown ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr); 159c4762a1bSJed Brown ierr = TSSetType(ts,TSBASICSYMPLECTIC);CHKERRQ(ierr); 160c4762a1bSJed Brown ierr = TSRHSSplitSetIS(ts,"position",is1);CHKERRQ(ierr); 161c4762a1bSJed Brown ierr = TSRHSSplitSetIS(ts,"momentum",is2);CHKERRQ(ierr); 162c4762a1bSJed Brown ierr = TSRHSSplitSetRHSFunction(ts,"position",NULL,RHSFunction1,&user);CHKERRQ(ierr); 163c4762a1bSJed Brown ierr = TSRHSSplitSetRHSFunction(ts,"momentum",NULL,RHSFunction2,&user);CHKERRQ(ierr); 164c4762a1bSJed Brown ierr = TSSetRHSFunction(ts,NULL,RHSFunction,&user);CHKERRQ(ierr); 165c4762a1bSJed Brown 166c4762a1bSJed Brown ierr = TSSetMaxTime(ts,ftime);CHKERRQ(ierr); 167c4762a1bSJed Brown ierr = TSSetTimeStep(ts,0.0001);CHKERRQ(ierr); 168c4762a1bSJed Brown ierr = TSSetMaxSteps(ts,1000);CHKERRQ(ierr); 169c4762a1bSJed Brown ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr); 170c4762a1bSJed Brown if (monitor) { 171c4762a1bSJed Brown ierr = TSMonitorSet(ts,Monitor,&user,NULL);CHKERRQ(ierr); 172c4762a1bSJed Brown } 173c4762a1bSJed Brown 174c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 175c4762a1bSJed Brown Set initial conditions 176c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 177c4762a1bSJed Brown ierr = VecGetArray(U,&u_ptr);CHKERRQ(ierr); 178c4762a1bSJed Brown u_ptr[0] = 0.2; 179c4762a1bSJed Brown u_ptr[1] = 0.0; 180c4762a1bSJed Brown ierr = VecRestoreArray(U,&u_ptr);CHKERRQ(ierr); 181c4762a1bSJed Brown 182c4762a1bSJed Brown ierr = TSSetTime(ts,0.0);CHKERRQ(ierr); 183c4762a1bSJed Brown 184c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 185c4762a1bSJed Brown Set runtime options 186c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 187c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 188c4762a1bSJed Brown 189c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 190c4762a1bSJed Brown Solve nonlinear system 191c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 192c4762a1bSJed Brown ierr = TSSolve(ts,U);CHKERRQ(ierr); 193c4762a1bSJed Brown ierr = TSGetSolveTime(ts,&ftime);CHKERRQ(ierr); 194c4762a1bSJed Brown ierr = VecView(U,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 195c4762a1bSJed Brown 196c4762a1bSJed Brown ierr = 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));CHKERRQ(ierr); 197c4762a1bSJed Brown 198c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 199c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 200c4762a1bSJed Brown are no longer needed. 201c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 202c4762a1bSJed Brown ierr = VecDestroy(&U);CHKERRQ(ierr); 203c4762a1bSJed Brown ierr = TSDestroy(&ts);CHKERRQ(ierr); 204c4762a1bSJed Brown ierr = ISDestroy(&is1);CHKERRQ(ierr); 205c4762a1bSJed Brown ierr = ISDestroy(&is2);CHKERRQ(ierr); 206c4762a1bSJed Brown ierr = PetscFinalize(); 207c4762a1bSJed Brown return ierr; 208c4762a1bSJed Brown } 209c4762a1bSJed Brown 210c4762a1bSJed Brown /*TEST 211c4762a1bSJed Brown build: 212c4762a1bSJed Brown requires: !single !complex 213c4762a1bSJed Brown 214c4762a1bSJed Brown test: 215c4762a1bSJed Brown args: -ts_basicsymplectic_type 1 -monitor 216c4762a1bSJed Brown 217c4762a1bSJed Brown test: 218c4762a1bSJed Brown suffix: 2 219c4762a1bSJed Brown args: -ts_basicsymplectic_type 2 -monitor 220c4762a1bSJed Brown 221c4762a1bSJed Brown TEST*/ 222