diff options
author | diener <diener@2a26948c-0e4f-0410-aee8-f1d3e353619c> | 2003-02-11 00:45:19 +0000 |
---|---|---|
committer | diener <diener@2a26948c-0e4f-0410-aee8-f1d3e353619c> | 2003-02-11 00:45:19 +0000 |
commit | d4514f75e0ac3be78f4dcfee1b0f23a540047142 (patch) | |
tree | e4b764105caf3860f72bead23aa17f956d40776e /src/EHFinder_FindSurface.F90 | |
parent | 3144b24a45f6a6279f83dddb787ee30336f101bc (diff) |
Fixed a bug that only showed up on 64 bit processors. An array of CCTK_Pointer
should have been an array of CCTK_INT. On my laptop it didn't make any
difference but on lemieux and seaborg it did.
git-svn-id: http://svn.einsteintoolkit.org/cactus/EinsteinAnalysis/EHFinder/trunk@92 2a26948c-0e4f-0410-aee8-f1d3e353619c
Diffstat (limited to 'src/EHFinder_FindSurface.F90')
-rw-r--r-- | src/EHFinder_FindSurface.F90 | 13 |
1 files changed, 6 insertions, 7 deletions
diff --git a/src/EHFinder_FindSurface.F90 b/src/EHFinder_FindSurface.F90 index d9be219..ece1b1c 100644 --- a/src/EHFinder_FindSurface.F90 +++ b/src/EHFinder_FindSurface.F90 @@ -31,10 +31,10 @@ subroutine EHFinder_FindSurface(CCTK_ARGUMENTS) crange_min_glob, crange_max_glob CCTK_REAL :: ltheta, utheta, lphi, uphi - CCTK_REAL :: dtheta, dphi, dthetainv, dphiinv, delta, rloc + CCTK_REAL :: dtheta, dphi, dthetainv, dphiinv, min_delta, rloc CCTK_POINTER, dimension(3) :: interp_coords CCTK_POINTER, dimension(4) :: out_array - CCTK_POINTER :: in_array + CCTK_INT :: in_array CCTK_INT, dimension(4), parameter :: op_indices = (/ 0, 0, 0, 0 /), & op_codes = (/ 0, 1, 2, 3 /), & out_types = (/ CCTK_VARIABLE_REAL, & @@ -93,7 +93,7 @@ subroutine EHFinder_FindSurface(CCTK_ARGUMENTS) ! IF not find an approximation for it by averaging the coodinates that ! are within one grid cell distance from the surface in question. nc_loc = 0; center_loc = zero - delta = minval ( cctk_delta_space ) + min_delta = minval ( cctk_delta_space ) do k = kzl, kzr do j = jyl, jyr do i = ixl, ixr @@ -363,7 +363,7 @@ subroutine EHFinder_FindSurface(CCTK_ARGUMENTS) ! a help in keeping track on how long ago drsurf for a given point has been ! reduced and is initialised to -1 to indicate that no reduction has taken ! place. - drsurf = two * delta + drsurf = two * min_delta rsurf = zero n_since_last_reduction = -1 @@ -373,8 +373,7 @@ subroutine EHFinder_FindSurface(CCTK_ARGUMENTS) ! call CCTK_InterpHandle ( interp_handle, "Lagrange polynomial interpolation" ) call CCTK_InterpHandle ( interp_handle, "Hermite polynomial interpolation" ) if ( interp_handle .lt. 0 ) then - call CCTK_WARN( 0, "Cannot get handle for interpolation. & - & Forgot to activate an implementation providing interpolation operators??" ) + call CCTK_WARN( 0, "Cannot get handle for interpolation. Forgot to activate an implementation providing interpolation operators??" ) end if ! For now order=2 is hard wired. @@ -446,7 +445,7 @@ subroutine EHFinder_FindSurface(CCTK_ARGUMENTS) ! print*,maxl,rsurf(maxl(1),maxl(2)),& ! drsurf(maxl(1),maxl(2)),& ! f_interp(maxl(1),maxl(2)) - if ( maxdr .lt. delta ) exit find_starting_points + if ( maxdr .lt. min_delta ) exit find_starting_points do j = 1, lsh(2) do i = 1, lsh(1) |