aboutsummaryrefslogtreecommitdiff
diff options
context:
space:
mode:
-rw-r--r--src/WaveToy.F24
1 files changed, 18 insertions, 6 deletions
diff --git a/src/WaveToy.F b/src/WaveToy.F
index 17022d3..27dea4a 100644
--- a/src/WaveToy.F
+++ b/src/WaveToy.F
@@ -37,7 +37,10 @@
! Declare local variables
INTEGER :: i,j,k
INTEGER :: istart, jstart, kstart, iend, jend, kend
+
CCTK_REAL :: dx,dy,dz,dt
+ CCTK_REAL :: dx2,dy2,dz2,dt2
+ CCTK_REAL :: dx2i,dy2i,dz2i
! Set up shorthands
! -----------------
@@ -46,6 +49,15 @@
dz = CCTK_DELTA_SPACE(3)
dt = CCTK_DELTA_TIME
+ dx2 = dx*dx
+ dy2 = dy*dy
+ dz2 = dz*dz
+ dt2 = dt*dt
+
+ dx2i = 1.0/dx2
+ dy2i = 1.0/dy2
+ dz2i = 1.0/dz2
+
istart = 2
jstart = 2
kstart = 2
@@ -61,12 +73,12 @@
do i = istart, iend
phi_tmp(i,j,k) =
- 1 2.0*(1.0 - (dt**2)*(1.0/dx**2 +
- 2 1.0/dy**2 +1.0/dz**2))*phi(i,j,k) -
- 3 phi_old(i,j,k) + (dt**2) *
- 5 ((phi(i+1,j,k)+phi(i-1,j,k))/dx**2
- 6 +(phi(i,j+1,k)+phi(i,j-1,k))/dy**2
- 7 +(phi(i,j,k+1)+phi(i,j,k-1))/dz**2)
+ 1 2.0*(1.0 - (dt2)*(dx2i +
+ 2 dy2i + dz2i))*phi(i,j,k) -
+ 3 phi_old(i,j,k) + (dt2) *
+ 5 ((phi(i+1,j,k)+phi(i-1,j,k))*dx2i
+ 6 +(phi(i,j+1,k)+phi(i,j-1,k))*dy2i
+ 7 +(phi(i,j,k+1)+phi(i,j,k-1))*dz2i)
end do
end do