aboutsummaryrefslogtreecommitdiff
path: root/src
diff options
context:
space:
mode:
authorknarf <knarf@c83d129a-5a75-4d5a-9c4d-ed3a5855bf45>2010-04-13 22:40:33 +0000
committerknarf <knarf@c83d129a-5a75-4d5a-9c4d-ed3a5855bf45>2010-04-13 22:40:33 +0000
commitdbe920a0139b02848ad74d8c139a48293be31efb (patch)
treeb521b861e0cc4d44667db37c8d150e5e768cd757 /src
parent5d8bb248ae7fa56ae89e404c1a04e5b8de62a6dd (diff)
correct logical/integer expressions
git-svn-id: http://svn.einsteintoolkit.org/cactus/EinsteinEvolve/GRHydro/trunk@106 c83d129a-5a75-4d5a-9c4d-ed3a5855bf45
Diffstat (limited to 'src')
-rw-r--r--src/Whisky_Con2Prim.F9016
-rw-r--r--src/Whisky_ENOReconstruct.F902
-rw-r--r--src/Whisky_PPM.F9012
-rw-r--r--src/Whisky_TVDReconstruct.F906
-rw-r--r--src/Whisky_UpdateMask.F902
5 files changed, 18 insertions, 20 deletions
diff --git a/src/Whisky_Con2Prim.F90 b/src/Whisky_Con2Prim.F90
index 344b138..f66eedb 100644
--- a/src/Whisky_Con2Prim.F90
+++ b/src/Whisky_Con2Prim.F90
@@ -604,7 +604,7 @@ subroutine Conservative2PrimitiveBounds(CCTK_ARGUMENTS)
!do not compute if in atmosphere or in an excised region
if (SpaceMask_CheckStateBitsF90(space_mask, i, j, k, type_bits, atmosphere) .or. &
- whisky_enable_internal_excision .and. (hydro_excision_mask(i,j,k) .ne. 0)) cycle
+ whisky_enable_internal_excision /= 0 .and. (hydro_excision_mask(i,j,k) .ne. 0)) cycle
gxxl = 0.5d0 * (gxx(i,j,k) + gxx(i-xoffset,j-yoffset,k-zoffset))
gxyl = 0.5d0 * (gxy(i,j,k) + gxy(i-xoffset,j-yoffset,k-zoffset))
@@ -818,7 +818,7 @@ subroutine Conservative2PrimitivePolytype(CCTK_ARGUMENTS)
!do not compute if in atmosphere or in an excised region
if (SpaceMask_CheckStateBitsF90(space_mask, i, j, k, type_bits, atmosphere) .or. &
- whisky_enable_internal_excision .and. (hydro_excision_mask(i,j,k) .ne. 0)) cycle
+ whisky_enable_internal_excision /= 0 .and. (hydro_excision_mask(i,j,k) .ne. 0)) cycle
if (conformal_state .eq. 0) then
call SpatialDeterminant(gxx(i,j,k),gxy(i,j,k),gxz(i,j,k),&
@@ -1146,7 +1146,7 @@ subroutine Con2PrimBoundsPolytype(CCTK_ARGUMENTS)
!do not compute if in atmosphere or in an excised region
if (SpaceMask_CheckStateBitsF90(space_mask, i, j, k, type_bits, atmosphere) .or. &
- whisky_enable_internal_excision .and. (hydro_excision_mask(i,j,k) .ne. 0)) cycle
+ whisky_enable_internal_excision /= 0 .and. (hydro_excision_mask(i,j,k) .ne. 0)) cycle
gxxl = 0.5d0 * (gxx(i,j,k) + gxx(i-xoffset,j-yoffset,k-zoffset))
gxyl = 0.5d0 * (gxy(i,j,k) + gxy(i-xoffset,j-yoffset,k-zoffset))
@@ -1629,7 +1629,7 @@ subroutine Conservative2PrimitiveGeneral(CCTK_ARGUMENTS)
!do not compute if in atmosphere or in an excised region
if (SpaceMask_CheckStateBitsF90(space_mask, i, j, k, type_bits, atmosphere) .or. &
- whisky_enable_internal_excision .and. (hydro_excision_mask(i,j,k) .ne. 0)) &
+ whisky_enable_internal_excision /= 0 .and. (hydro_excision_mask(i,j,k) .ne. 0)) &
then
press_new(i,j,k) = press(i,j,k)
press_old(i,j,k) = press(i,j,k)
@@ -1733,7 +1733,7 @@ subroutine Conservative2PrimitiveGeneral(CCTK_ARGUMENTS)
!do not compute if in atmosphere or in an excised region
if (SpaceMask_CheckStateBitsF90(space_mask, i, j, k, type_bits, atmosphere) .or. &
- whisky_enable_internal_excision .and. (hydro_excision_mask(i,j,k) .ne. 0)) cycle
+ whisky_enable_internal_excision /= 0 .and. (hydro_excision_mask(i,j,k) .ne. 0)) cycle
f = press_new(i,j,k) - press(i,j,k)
@@ -1787,7 +1787,7 @@ subroutine Conservative2PrimitiveGeneral(CCTK_ARGUMENTS)
!do not compute if in atmosphere or in an excised region
if (SpaceMask_CheckStateBitsF90(space_mask, i, j, k, type_bits, atmosphere) .or. &
- whisky_enable_internal_excision .and. (hydro_excision_mask(i,j,k) .ne. 0)) cycle
+ whisky_enable_internal_excision /= 0 .and. (hydro_excision_mask(i,j,k) .ne. 0)) cycle
!!$ Note that in the following we enforce eps > 0
!!$ This makes the later check redundant
@@ -2017,7 +2017,7 @@ subroutine Con2PrimPolytypeGeneral(CCTK_ARGUMENTS)
!do not compute if in atmosphere or in an excised region
if (SpaceMask_CheckStateBitsF90(space_mask, i, j, k, type_bits, atmosphere) .or. &
- whisky_enable_internal_excision .and. (hydro_excision_mask(i,j,k) .ne. 0)) cycle
+ whisky_enable_internal_excision /= 0 .and. (hydro_excision_mask(i,j,k) .ne. 0)) cycle
enthalpy = 1.d0 + eps(i,j,k) + press(i,j,k) / press_new(i,j,k)
@@ -2070,7 +2070,7 @@ subroutine Con2PrimPolytypeGeneral(CCTK_ARGUMENTS)
!do not compute if in atmosphere or in an excised region
if (SpaceMask_CheckStateBitsF90(space_mask, i, j, k, type_bits, atmosphere) .or. &
- whisky_enable_internal_excision .and. (hydro_excision_mask(i,j,k) .ne. 0)) cycle
+ whisky_enable_internal_excision /= 0 .and. (hydro_excision_mask(i,j,k) .ne. 0)) cycle
udens = dens(i,j,k) / sqrt(Whisky_det(i,j,k))
diff --git a/src/Whisky_ENOReconstruct.F90 b/src/Whisky_ENOReconstruct.F90
index 4fc699c..6e73b21 100644
--- a/src/Whisky_ENOReconstruct.F90
+++ b/src/Whisky_ENOReconstruct.F90
@@ -181,7 +181,7 @@ subroutine Whisky_ENOReconstruct1d(order, nx, v, vminus, vplus, trivial_rp, &
!!$ Initialize excision
do i = 1, nx
- if (whisky_enable_internal_excision .and. (hydro_excision_mask(i) .ne. 0)) then
+ if (whisky_enable_internal_excision /= 0 .and. (hydro_excision_mask(i) .ne. 0)) then
vdiff(1, i) = large * (-1.d0*(1.d0+mod(i,10)))**i
trivial_rp(i) = .true.
excise(i) = .true.
diff --git a/src/Whisky_PPM.F90 b/src/Whisky_PPM.F90
index f04e535..4dbeb97 100644
--- a/src/Whisky_PPM.F90
+++ b/src/Whisky_PPM.F90
@@ -430,9 +430,8 @@ do i = whisky_stencil, nx - whisky_stencil + 1
end if
!!$ excision
- if (hydro_excision) then
do i = 1, nx
- if (whisky_enable_internal_excision .and. &
+ if (whisky_enable_internal_excision /= 0 .and. &
(hydro_excision_mask(i) .ne. 0)) then
if (i .gt. 1) then
trivial_rp(i-1)=.true.
@@ -443,7 +442,7 @@ do i = whisky_stencil, nx - whisky_stencil + 1
!!$ have to follow the order of sub-expressions given here and might
!!$ access outside the array range
cond = i .gt. 1
- if (cond .and. whisky_enable_internal_excision) then
+ if (cond .and. whisky_enable_internal_excision /= 0) then
cond = hydro_excision_mask(i-1) .ne. 0
end if
if (cond) then
@@ -471,7 +470,7 @@ do i = whisky_stencil, nx - whisky_stencil + 1
end if
else
cond = (i.gt.2) .and. (i.lt.nx)
- if (cond .and. whisky_enable_internal_excision) then
+ if (cond .and. whisky_enable_internal_excision /= 0) then
cond = (ppm_mppm .eq. 0) .and. (hydro_excision_mask(i-2) .ne. 0)
end if
if (cond) then
@@ -485,7 +484,7 @@ do i = whisky_stencil, nx - whisky_stencil + 1
end if
end if
cond = i.lt.nx
- if (cond .and. whisky_enable_internal_excision) then
+ if (cond .and. whisky_enable_internal_excision /= 0) then
cond = hydro_excision_mask(i+1) .ne. 0
end if
if (cond) then
@@ -513,7 +512,7 @@ do i = whisky_stencil, nx - whisky_stencil + 1
endif
else
cond = (i.lt.nx-1) .and. (i.gt.1)
- if (cond .and. whisky_enable_internal_excision) then
+ if (cond .and. whisky_enable_internal_excision /= 0) then
cond = (ppm_mppm .eq. 0) .and. (hydro_excision_mask(i+2) .ne. 0)
end if
if (cond) then
@@ -528,7 +527,6 @@ do i = whisky_stencil, nx - whisky_stencil + 1
end if
end if
end do
- end if
return
end subroutine SimplePPM_1d
diff --git a/src/Whisky_TVDReconstruct.F90 b/src/Whisky_TVDReconstruct.F90
index 3141fb4..4459da2 100644
--- a/src/Whisky_TVDReconstruct.F90
+++ b/src/Whisky_TVDReconstruct.F90
@@ -56,18 +56,18 @@ subroutine tvdreconstruct(nx, ny, nz, xoffset, yoffset, zoffset, &
do k = whisky_stencil, nz-whisky_stencil+1
do j = whisky_stencil, ny-whisky_stencil+1
do i = whisky_stencil, nx-whisky_stencil+1
- if (whisky_enable_internal_excision .and. &
+ if (whisky_enable_internal_excision /= 0 .and. &
(hydro_excision_mask(i,j,k) .ne. 0)) then
trivial_rp(i-xoffset, j-yoffset, k-zoffset) = .true.
trivial_rp(i, j, k) = .true.
bextm(i, j, k) = orig(i, j, k)
bextp(i, j, k) = orig(i, j, k)
- if (whisky_enable_internal_excision .and. &
+ if (whisky_enable_internal_excision /= 0 .and. &
(hydro_excision_mask(i+xoffset,j+yoffset,k+zoffset) .eq. 0)) then
bextm(i, j, k) = orig(i+xoffset, j+yoffset, k+zoffset)
bextp(i, j, k) = orig(i+xoffset, j+yoffset, k+zoffset)
end if
- else if (whisky_enable_internal_excision .and. &
+ else if (whisky_enable_internal_excision /= 0 .and. &
((hydro_excision_mask(i-xoffset,j-yoffset,k-zoffset) .ne. 0) .or. &
(hydro_excision_mask(i+xoffset,j+yoffset,k+zoffset) .ne. 0))) then
bextm(i, j, k) = orig(i, j, k)
diff --git a/src/Whisky_UpdateMask.F90 b/src/Whisky_UpdateMask.F90
index 5d79424..d67600f 100644
--- a/src/Whisky_UpdateMask.F90
+++ b/src/Whisky_UpdateMask.F90
@@ -38,7 +38,7 @@ subroutine WhiskyUpdateAtmosphereMask(CCTK_ARGUMENTS)
do k = 1, cctk_lsh(3)
do j = 1, cctk_lsh(2)
do i = 1, cctk_lsh(1)
- if ( whisky_enable_internal_excision .and. (hydro_excision_mask(i,j,k) .ne. 0) .or. &
+ if ( whisky_enable_internal_excision /= 0 .and. (hydro_excision_mask(i,j,k) .ne. 0) .or. &
(SpaceMask_CheckStateBitsF90(space_mask, i, j, k, \
type_bits, atmosphere)) .or. &
(tau(i,j,k) + frac * taurhs(i,j,k) .le. 0.d0) .or. &