From e07ae62012ec2ccb3d283670a90252ea6d4bf6cd Mon Sep 17 00:00:00 2001 From: goodale Date: Tue, 28 Mar 2000 07:04:08 +0000 Subject: 60% more MFlop/s Tom git-svn-id: http://svn.cactuscode.org/arrangements/CactusWave/WaveToyFreeF90/trunk@15 daab82bb-f315-4ad1-b6d0-9353ff8b6e27 --- src/WaveToy.F90 | 24 ++++++++++++++++++------ 1 file changed, 18 insertions(+), 6 deletions(-) diff --git a/src/WaveToy.F90 b/src/WaveToy.F90 index eaa279c..1a02912 100644 --- a/src/WaveToy.F90 +++ b/src/WaveToy.F90 @@ -42,7 +42,10 @@ subroutine WaveToyFreeF90_Evolution(CCTK_ARGUMENTS) ! 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 ! ----------------- @@ -51,6 +54,15 @@ subroutine WaveToyFreeF90_Evolution(CCTK_ARGUMENTS) 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 @@ -66,12 +78,12 @@ subroutine WaveToyFreeF90_Evolution(CCTK_ARGUMENTS) do i = istart, iend phi_tmp(i,j,k) = & - 2.0*(1.0 - (dt**2)*(1.0/dx**2 + & - 1.0/dy**2 +1.0/dz**2))*phi(i,j,k) - & - phi_old(i,j,k) + (dt**2) * & - ((phi(i+1,j,k)+phi(i-1,j,k))/dx**2 & - +(phi(i,j+1,k)+phi(i,j-1,k))/dy**2 & - +(phi(i,j,k+1)+phi(i,j,k-1))/dz**2) + 2.0*(1.0 - (dt2)*(dx2i + & + dy2i + dz2i))*phi(i,j,k) - & + phi_old(i,j,k) + (dt2) * & + ((phi(i+1,j,k)+phi(i-1,j,k))*dx2i & + +(phi(i,j+1,k)+phi(i,j-1,k))*dy2i & + +(phi(i,j,k+1)+phi(i,j,k-1))*dz2i) end do end do -- cgit v1.2.3