aboutsummaryrefslogtreecommitdiff
path: root/src
diff options
context:
space:
mode:
authordiener <diener@2a26948c-0e4f-0410-aee8-f1d3e353619c>2003-04-15 07:24:37 +0000
committerdiener <diener@2a26948c-0e4f-0410-aee8-f1d3e353619c>2003-04-15 07:24:37 +0000
commit4aafcc75d35b54d5e4d23cb7b7da771023ba80f6 (patch)
treefe8b9af9450cfb6868f2b7d162bf0ea614bc5a7a /src
parentde02c721377d8c76059d87c9dfa8676499dfb272 (diff)
Fixed a bug that made the surface finding not work properly on
multiple processors. git-svn-id: http://svn.einsteintoolkit.org/cactus/EinsteinAnalysis/EHFinder/trunk@95 2a26948c-0e4f-0410-aee8-f1d3e353619c
Diffstat (limited to 'src')
-rw-r--r--src/EHFinder_FindSurface.F90206
1 files changed, 110 insertions, 96 deletions
diff --git a/src/EHFinder_FindSurface.F90 b/src/EHFinder_FindSurface.F90
index 655109b..8f1f1e7 100644
--- a/src/EHFinder_FindSurface.F90
+++ b/src/EHFinder_FindSurface.F90
@@ -161,113 +161,126 @@ subroutine EHFinder_FindSurface(CCTK_ARGUMENTS)
sym_factor = one
s_symx = .false.; s_symy = .false.; s_symz = .false.
- ! Bitant with symmetry in the z-direction
- if ( ( .not. symx ) .and. ( .not. symy ) .and. ( symz ) ) then
- if ( crange_min_glob(3) .lt. zero ) then
- ltheta = zero; utheta = half * pi
- lphi = zero; uphi = two * pi
- sym_factor = two
- center(3) = zero
- s_symz = .true.
- end if
- end if
+ if ( CCTK_EQUALS ( domain, 'bitant' ) ) then
+
+ ! Bitant with symmetry in the z-direction
+
+ if ( CCTK_EQUALS ( bitant_plane, 'xy' ) ) then
+ if ( crange_min_glob(3) .lt. zero ) then
+ ltheta = zero; utheta = half * pi
+ lphi = zero; uphi = two * pi
+ sym_factor = two
+ center(3) = zero
+ s_symz = .true.
+ end if
+ end if
- ! Bitant with symmetry in the y-direction
- if ( ( .not. symx ) .and. ( symy ) .and. ( .not. symz ) ) then
- if ( crange_min_glob(2) .lt. zero ) then
- ltheta = zero; utheta = pi
- lphi = zero; uphi = pi
- sym_factor = two
- center(2) = zero
- s_symy = .true.
- end if
- end if
+ ! Bitant with symmetry in the y-direction
+
+ if ( CCTK_EQUALS ( bitant_plane, 'xz' ) ) then
+ if ( crange_min_glob(2) .lt. zero ) then
+ ltheta = zero; utheta = pi
+ lphi = zero; uphi = pi
+ sym_factor = two
+ center(2) = zero
+ s_symy = .true.
+ end if
+ end if
- ! Bitant with symmetry in the x-direction
- if ( ( symx ) .and. ( .not. symy ) .and. ( .not. symz ) ) then
- if ( crange_min_glob(1) .lt. zero ) then
- ltheta = zero; utheta = pi
- lphi = -half * pi; uphi = half * pi
- sym_factor = two
- center(1) = zero
- s_symx = .true.
- end if
- end if
+ ! Bitant with symmetry in the x-direction
- ! Quadrant in x-direction
- if ( ( .not. symx ) .and. ( symy ) .and. ( symz ) ) then
- if ( ( crange_min_glob(3) .lt. zero ) .and. &
- ( crange_min_glob(2) .lt. zero ) ) then
- ltheta = zero; utheta = half * pi
- lphi = zero; uphi = pi
- sym_factor = four
- center(3) = zero; center(2) = zero
- s_symz = .true.; s_symy = .true.
- else if ( crange_min_glob(3) .lt. zero ) then
- ltheta = zero; utheta = half * pi
- lphi = zero; uphi = two * pi
- sym_factor = two
- center(3) = zero
- s_symz = .true.
- else if ( crange_min_glob(2) .lt. zero ) then
- ltheta = zero; utheta = pi
- lphi = zero; uphi = pi
- sym_factor = two
- center(2) = zero
- s_symy = .true.
+ if ( CCTK_EQUALS ( bitant_plane, 'yz' ) ) then
+ if ( crange_min_glob(1) .lt. zero ) then
+ ltheta = zero; utheta = pi
+ lphi = -half * pi; uphi = half * pi
+ sym_factor = two
+ center(1) = zero
+ s_symx = .true.
+ end if
end if
end if
+
+ if ( CCTK_EQUALS ( domain, 'quadrant' ) ) then
+
+ ! Quadrant in x-direction
+
+ if ( CCTK_EQUALS ( quadrant_direction, 'x' ) ) then
+ if ( ( crange_min_glob(3) .lt. zero ) .and. &
+ ( crange_min_glob(2) .lt. zero ) ) then
+ ltheta = zero; utheta = half * pi
+ lphi = zero; uphi = pi
+ sym_factor = four
+ center(3) = zero; center(2) = zero
+ s_symz = .true.; s_symy = .true.
+ else if ( crange_min_glob(3) .lt. zero ) then
+ ltheta = zero; utheta = half * pi
+ lphi = zero; uphi = two * pi
+ sym_factor = two
+ center(3) = zero
+ s_symz = .true.
+ else if ( crange_min_glob(2) .lt. zero ) then
+ ltheta = zero; utheta = pi
+ lphi = zero; uphi = pi
+ sym_factor = two
+ center(2) = zero
+ s_symy = .true.
+ end if
+ end if
- ! Quadrant in y-direction
- if ( ( symx ) .and. ( .not. symy ) .and. ( symz ) ) then
- if ( ( crange_min_glob(3) .lt. zero ) .and. &
- ( crange_min_glob(1) .lt. zero ) ) then
- ltheta = zero; utheta = half * pi
- lphi = -half * pi; uphi = half * pi
- sym_factor = four
- center(3) = zero; center(1) = zero
- s_symz = .true.; s_symx = .true.
- else if ( crange_min_glob(3) .lt. zero ) then
- ltheta = zero; utheta = half * pi
- lphi = zero; uphi = two * pi
- sym_factor = two
- center(3) = zero
- s_symz = .true.
- else if ( crange_min_glob(1) .lt. zero ) then
- ltheta = zero; utheta = pi
- lphi = -half * pi; uphi = half * pi
- sym_factor = two
- center(1) = zero
- s_symx = .true.
+ ! Quadrant in y-direction
+
+ if ( CCTK_EQUALS ( quadrant_direction, 'y' ) ) then
+ if ( ( crange_min_glob(3) .lt. zero ) .and. &
+ ( crange_min_glob(1) .lt. zero ) ) then
+ ltheta = zero; utheta = half * pi
+ lphi = -half * pi; uphi = half * pi
+ sym_factor = four
+ center(3) = zero; center(1) = zero
+ s_symz = .true.; s_symx = .true.
+ else if ( crange_min_glob(3) .lt. zero ) then
+ ltheta = zero; utheta = half * pi
+ lphi = zero; uphi = two * pi
+ sym_factor = two
+ center(3) = zero
+ s_symz = .true.
+ else if ( crange_min_glob(1) .lt. zero ) then
+ ltheta = zero; utheta = pi
+ lphi = -half * pi; uphi = half * pi
+ sym_factor = two
+ center(1) = zero
+ s_symx = .true.
+ end if
end if
- end if
- ! Quadrant in z-direction
- if ( ( symx ) .and. ( symy ) .and. ( .not. symz ) ) then
- if ( ( crange_min_glob(2) .lt. zero ) .and. &
- ( crange_min_glob(1) .lt. zero ) ) then
- ltheta = zero; utheta = pi
- lphi = zero; uphi = half * pi
- sym_factor = four
- center(2) = zero; center(1) = zero
- s_symy = .true.; s_symx = .true.
- else if ( crange_min_glob(2) .lt. zero ) then
- ltheta = zero; utheta = pi
- lphi = zero; uphi = pi
- sym_factor = two
- center(2) = zero
- s_symy = .true.
- else if ( crange_min_glob(1) .lt. zero ) then
- ltheta = zero; utheta = pi
- lphi = -half * pi; uphi = half * pi
- sym_factor = two
- center(1) = zero
- s_symx = .true.
+ ! Quadrant in z-direction
+
+ if ( CCTK_EQUALS ( quadrant_direction,'z' ) ) then
+ if ( ( crange_min_glob(2) .lt. zero ) .and. &
+ ( crange_min_glob(1) .lt. zero ) ) then
+ ltheta = zero; utheta = pi
+ lphi = zero; uphi = half * pi
+ sym_factor = four
+ center(2) = zero; center(1) = zero
+ s_symy = .true.; s_symx = .true.
+ else if ( crange_min_glob(2) .lt. zero ) then
+ ltheta = zero; utheta = pi
+ lphi = zero; uphi = pi
+ sym_factor = two
+ center(2) = zero
+ s_symy = .true.
+ else if ( crange_min_glob(1) .lt. zero ) then
+ ltheta = zero; utheta = pi
+ lphi = -half * pi; uphi = half * pi
+ sym_factor = two
+ center(1) = zero
+ s_symx = .true.
+ end if
end if
end if
! Octant
- if ( ( symx ) .and. ( symy ) .and. ( symz ) ) then
+
+ if ( CCTK_EQUALS ( domain, 'octant' ) ) then
if ( ( crange_min_glob(3) .lt. zero ) .and. &
( crange_min_glob(2) .lt. zero ) .and. &
( crange_min_glob(1) .lt. zero ) ) then
@@ -337,6 +350,7 @@ subroutine EHFinder_FindSurface(CCTK_ARGUMENTS)
! print*,ltheta/pi,utheta/pi,lphi/pi,uphi/pi
! print*,sym_factor
! print*,center
+! print*,symx,symy,symz
! ! For now symmetries are not handled properly
! if ( symx ) center(1) = zero
! if ( symy ) center(2) = zero
@@ -648,7 +662,7 @@ subroutine EHFinder_FindSurface(CCTK_ARGUMENTS)
! print*,'Calculated derivatives'
! print*,rsurf
call Util_TableDestroy ( status, table_handle )
- print*,'gjphi = ',gjphi
+! print*,'gjphi = ',gjphi
end subroutine EHFinder_FindSurface