aboutsummaryrefslogtreecommitdiff
diff options
context:
space:
mode:
authorpollney <pollney@89daf98e-ef62-4674-b946-b8ff9de2216c>2003-01-06 15:55:44 +0000
committerpollney <pollney@89daf98e-ef62-4674-b946-b8ff9de2216c>2003-01-06 15:55:44 +0000
commitb9c951cd037aabc734f625244e60a5311ade7f87 (patch)
treec1f74afad335e9639580ae2635c58ea37e2a7f1f
parentc58747e1728af63ac429e3f97e7107f8fb6828d1 (diff)
Modified ahf.gauss and *.alm headers so that output format is identical with the test suites.
git-svn-id: http://svn.einsteintoolkit.org/cactus/EinsteinAnalysis/AHFinder/trunk@323 89daf98e-ef62-4674-b946-b8ff9de2216c
-rw-r--r--src/AHFinder.F96
-rw-r--r--src/AHFinder_output.F28
-rw-r--r--src/SetSym.F118
3 files changed, 130 insertions, 112 deletions
diff --git a/src/AHFinder.F b/src/AHFinder.F
index 0d646a6..8b25e30 100644
--- a/src/AHFinder.F
+++ b/src/AHFinder.F
@@ -444,98 +444,8 @@
nonaxi = .false.
end if
-! Find out reflection symmetries depending on
-! the values of the parameters.
- if ((CCTK_Equals(ahf_octant,"yes").eq.1).or.
- . (CCTK_Equals(ahf_octant,"high").eq.1)) then
-
- refx = .true.
- refy = .true.
- refz = .true.
-
- else
-
- if (ahf_refx.ne.0) then
- refx = .true.
- else
- refx = .false.
- end if
-
- if (ahf_refy.ne.0) then
- refy = .true.
- else
- refy = .false.
- end if
-
- if (ahf_refz.ne.0) then
- refz = .true.
- else
- refz = .false.
- end if
-
- end if
-
- if (ahf_phi.eq.0) then
- refx = .true.
- refy = .true.
- end if
-
-! Force grid symmetries. If the grid has some
-! symmetry (octant, quadrant, bitant), and the
-! horizon is centered on one of the symmetry
-! planes, then the grid symmetries are forced
-! on it regardless of the values of the symmetry
-! parameters.
-
- if (CCTK_Equals(domain,"octant").eq.1) then
-
- if (xc.eq.zero) refx = .true.
- if (yc.eq.zero) refy = .true.
- if (zc.eq.zero) refz = .true.
-
- else if (CCTK_Equals(domain,"quadrant").eq.1) then
-
- if (xc.eq.zero) refx = .true.
- if (yc.eq.zero) refy = .true.
-
- else if (CCTK_Equals(domain,"bitant").eq.1) then
-
- if (CCTK_Equals(bitant_plane,"xy").eq.1) then
- if (zc.eq.zero) refz = .true.
- else if (CCTK_Equals(bitant_plane,"xz").eq.1) then
- if (yc.eq.zero) refy = .true.
- else
- if (xc.eq.zero) refx = .true.
- end if
-
- end if
-
-! Check if we are using cartoon. If we are, then the
-! corresponding symmetries are forced on the horizon.
-
- if (ahf_cartoon.ne.0) then
- cartoon = .true.
- else
- cartoon = .false.
- end if
-
- if (cartoon) then
- nonaxi = .false.
- refx = .true.
- refy = .true.
- end if
-
-! If desired, write values of symmetry flags to screen.
-
- if (veryver) then
- write(*,*)
- write(*,*) 'Symmetries used for horizon:'
- write(*,*) 'refx = ',refx
- write(*,*) 'refy = ',refy
- write(*,*) 'refz = ',refz
- write(*,*)
- end if
+ call AHFinder_SetReflections(CCTK_ARGUMENTS)
! Find the values of {stepx,stepy,stepz} depending
! on the symmetries.
@@ -1620,7 +1530,3 @@
! ***************
end subroutine AHFinder
-
-
-
-
diff --git a/src/AHFinder_output.F b/src/AHFinder_output.F
index c9f4f8b..9073c3b 100644
--- a/src/AHFinder_output.F
+++ b/src/AHFinder_output.F
@@ -131,6 +131,8 @@ c some compilers cannot trim an empty string, so we add at least one char
character(len=200) :: filename
+ call AHFinder_SetReflections(CCTK_ARGUMENTS)
+
open(1,file=filename,form='formatted', status='replace')
write(1,"(A20)") '# GAUSSIAN CURVATURE'
@@ -158,13 +160,8 @@ c some compilers cannot trim an empty string, so we add at least one char
write(1,"(A9,L1)") '# refy = ',refy
write(1,"(A9,L1)") '# refz = ',refz
write(1,"(A1)") '#'
- write(1,"(A11,I7)") '# ntheta = ',(ntheta+1)
- write(1,"(A11,I7)") '# nphi = ',(nphi+1)
-
- write(1,*)
- write(1,"(A11,I4)") '# Time step',cctk_iteration
- write(1,"(A6,ES11.3)") '# Time',cctk_time
- write(1,"(A6,I4)") '# Call',ahf_ncall
+ write(1,"(A11,I7)") '# ntheta = ',(ahf_ntheta+1)
+ write(1,"(A11,I7)") '# nphi = ',(ahf_nphi+1)
close(1)
@@ -458,21 +455,12 @@ c some compilers cannot trim an empty string, so we add at least one char
inquire(file=almf, exist=almf_exists)
if (almf_exists) then
open(1,file=almf,form='formatted',status='old',position='append')
- write(1,*)
else
open(1,file=almf,form='formatted',status='replace')
write(1,"(A21)") '# Radial coefficients'
write(1,"(A1)") '#'
end if
- if (status.and.report) then
- out_centerx = xc
- out_centery = yc
- out_centerz = zc
-
- write(1,"(A15,3ES14.6)") '# centered on: ', out_centerx, out_centery, out_centerz
- endif
-
write(1,"(A12,I4)") '# Time step ',cctk_iteration
write(1,"(A7,ES14.6)") '# Time ',cctk_time
write(1,"(A7,I4)") '# Call ',ahf_ncall
@@ -503,6 +491,12 @@ c some compilers cannot trim an empty string, so we add at least one char
write(1,"(A30,I4)") '# Surface found: Not a horizon'
end if
end if
+
+ out_centerx = xc
+ out_centery = yc
+ out_centerz = zc
+
+ write(1,"(A15,3ES14.6)") '# centered on: ', out_centerx, out_centery, out_centerz
write(1,"(A1)") '#'
write(1,"(A22)") '# a_lm l m'
@@ -525,7 +519,7 @@ c some compilers cannot trim an empty string, so we add at least one char
end do
end do
end if
-
+ write(1,*)
else
write(1,"(A18,I4)") '# No surface found'
end if
diff --git a/src/SetSym.F b/src/SetSym.F
index a7d4f45..a7c331b 100644
--- a/src/SetSym.F
+++ b/src/SetSym.F
@@ -77,3 +77,121 @@ c +,+,-
c End.
end subroutine AHFinder_SetSym
+
+
+ /*@@
+ @file SetSym.F
+ @date Jan 2003
+ @author Denis Pollney
+ @desc
+ Sets the reflection symmetry flags.
+ @enddesc
+ @@*/
+
+ subroutine AHFinder_SetReflections(CCTK_ARGUMENTS)
+ use AHFinder_dat
+
+ implicit none
+
+ DECLARE_CCTK_ARGUMENTS
+ DECLARE_CCTK_PARAMETERS
+ DECLARE_CCTK_FUNCTIONS
+
+ CCTK_REAL, parameter :: zero = 0.0D0
+ CCTK_REAL, parameter :: half = 0.5D0
+ CCTK_REAL, parameter :: one = 1.0D0
+
+! Find out reflection symmetries depending on
+! the values of the parameters.
+
+ if ((CCTK_Equals(ahf_octant,"yes").eq.1).or.
+ . (CCTK_Equals(ahf_octant,"high").eq.1)) then
+
+ refx = .true.
+ refy = .true.
+ refz = .true.
+
+ else
+
+ if (ahf_refx.ne.0) then
+ refx = .true.
+ else
+ refx = .false.
+ end if
+
+ if (ahf_refy.ne.0) then
+ refy = .true.
+ else
+ refy = .false.
+ end if
+
+ if (ahf_refz.ne.0) then
+ refz = .true.
+ else
+ refz = .false.
+ end if
+
+ end if
+
+ if (ahf_phi.eq.0) then
+ refx = .true.
+ refy = .true.
+ end if
+
+! Force grid symmetries. If the grid has some
+! symmetry (octant, quadrant, bitant), and the
+! horizon is centered on one of the symmetry
+! planes, then the grid symmetries are forced
+! on it regardless of the values of the symmetry
+! parameters.
+
+ if (CCTK_Equals(domain,"octant").eq.1) then
+
+ if (xc.eq.zero) refx = .true.
+ if (yc.eq.zero) refy = .true.
+ if (zc.eq.zero) refz = .true.
+
+ else if (CCTK_Equals(domain,"quadrant").eq.1) then
+
+ if (xc.eq.zero) refx = .true.
+ if (yc.eq.zero) refy = .true.
+
+ else if (CCTK_Equals(domain,"bitant").eq.1) then
+
+ if (CCTK_Equals(bitant_plane,"xy").eq.1) then
+ if (zc.eq.zero) refz = .true.
+ else if (CCTK_Equals(bitant_plane,"xz").eq.1) then
+ if (yc.eq.zero) refy = .true.
+ else
+ if (xc.eq.zero) refx = .true.
+ end if
+
+ end if
+
+! Check if we are using cartoon. If we are, then the
+! corresponding symmetries are forced on the horizon.
+
+ if (ahf_cartoon.ne.0) then
+ cartoon = .true.
+ else
+ cartoon = .false.
+ end if
+
+ if (cartoon) then
+ nonaxi = .false.
+ refx = .true.
+ refy = .true.
+ end if
+
+! If desired, write values of symmetry flags to screen.
+
+ if (veryver) then
+ write(*,*)
+ write(*,*) 'Symmetries used for horizon:'
+ write(*,*) 'refx = ',refx
+ write(*,*) 'refy = ',refy
+ write(*,*) 'refz = ',refz
+ write(*,*)
+ end if
+
+ end subroutine AHFinder_SetReflections