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