aboutsummaryrefslogtreecommitdiff
diff options
context:
space:
mode:
authorrhaas <rhaas@ef6f4158-a8ee-47d1-ba14-cb73256398e0>2012-01-25 15:19:00 +0000
committerrhaas <rhaas@ef6f4158-a8ee-47d1-ba14-cb73256398e0>2012-01-25 15:19:00 +0000
commitcc396ab996df8be585ddb70c84c9b943002c99a5 (patch)
treea1993a6e071038f155cd6623a541723b2e0c4e01
parentec2893918d69b8a3e5d9680a85fb540acfb49942 (diff)
support named spherical surfaces
git-svn-id: https://svn.cct.lsu.edu/repos/numrel/LSUThorns/QuasiLocalMeasures/trunk@36 ef6f4158-a8ee-47d1-ba14-cb73256398e0
-rw-r--r--interface.ccl5
-rw-r--r--param.ccl6
-rw-r--r--src/qlm_calculate.F902
-rw-r--r--src/qlm_import_surface.F9011
-rw-r--r--src/qlm_paramcheck.F907
-rw-r--r--test/qlm-minkowski.par3
6 files changed, 24 insertions, 10 deletions
diff --git a/interface.ccl b/interface.ccl
index 3cb4cd0..4510fce 100644
--- a/interface.ccl
+++ b/interface.ccl
@@ -10,6 +10,11 @@ CCTK_POINTER_TO_CONST \
FUNCTION GetMPICommWorld (CCTK_POINTER_TO_CONST IN cctkGH)
USES FUNCTION GetMPICommWorld
+# translate SphericalSurface name into ID
+CCTK_INT \
+FUNCTION sf_IdFromName (CCTK_INT IN id, CCTK_POINTER_TO_CONST IN name)
+USES FUNCTION sf_IdFromName
+
INT qlm_state[num_surfaces] TYPE=scalar
diff --git a/param.ccl b/param.ccl
index e40e623..ea78903 100644
--- a/param.ccl
+++ b/param.ccl
@@ -23,6 +23,12 @@ INT surface_index[100] "Spherical surface that contains the surface shape"
0:* :: "surface index"
} -1
+STRING surface_name[100] "Spherical surface that contains the surface shape"
+{
+ "" :: "use surface_index"
+ ".*" :: "surface name"
+} ""
+
STRING coordsystem "The coordinate system to use" STEERABLE=always
{
"" :: "must be a registered coordinate system"
diff --git a/src/qlm_calculate.F90 b/src/qlm_calculate.F90
index e9b602d..44574c6 100644
--- a/src/qlm_calculate.F90
+++ b/src/qlm_calculate.F90
@@ -55,7 +55,7 @@ subroutine qlm_calculate (CCTK_ARGUMENTS)
end if
if (hn > 0) then
- if (surface_index(hn) == -1) then
+ if (surface_index(hn) == -1 .and. CCTK_EQUALS(surface_name(hn), "")) then
qlm_calc_error(hn) = 1
qlm_have_valid_data(hn) = 0
qlm_have_killing_vector(hn) = 0
diff --git a/src/qlm_import_surface.F90 b/src/qlm_import_surface.F90
index d10b736..18173db 100644
--- a/src/qlm_import_surface.F90
+++ b/src/qlm_import_surface.F90
@@ -21,18 +21,19 @@ subroutine qlm_import_surface (CCTK_ARGUMENTS, hn)
call CCTK_INFO ("Importing surface shape")
end if
- if (surface_index(hn) < 0 .or. surface_index(hn) >= nsurfaces) then
+ if ((surface_index(hn) < 0 .or. surface_index(hn) >= nsurfaces) .and. &
+ CCTK_EQUALS(surface_name(hn), "")) then
call CCTK_WARN (0, "Illegal spherical surface index specified")
end if
+ sn = sf_IdFromName(surface_index(hn), surface_name(hn)) + 1
+
if (verbose/=0 .or. veryverbose/=0) then
- write (msg, '("Importing from spherical surface ",i4)') surface_index(hn)
+ write (msg, '("Importing from spherical surface ",a," ",i4)') &
+ surface_name(hn), sn - 1
call CCTK_INFO (msg)
end if
- sn = surface_index(hn) + 1
-
-
if (qlm_calc_error(hn) == 0 .and. cctk_iteration > qlm_iteration(hn)) then
diff --git a/src/qlm_paramcheck.F90 b/src/qlm_paramcheck.F90
index b8a8b8e..2e910a8 100644
--- a/src/qlm_paramcheck.F90
+++ b/src/qlm_paramcheck.F90
@@ -22,18 +22,19 @@ subroutine qlm_paramcheck (CCTK_ARGUMENTS)
do hn = 1, num_surfaces
- if (surface_index(hn) == -1) then
+ if (surface_index(hn) == -1 .and. CCTK_EQUALS(surface_name(hn),"")) then
! no surface selected; everything is fine
goto 9999
end if
- if (surface_index(hn) < 0 .or. surface_index(hn) >= nsurfaces) then
+ if ((surface_index(hn) < 0 .or. surface_index(hn) >= nsurfaces) .and. &
+ CCTK_EQUALS(surface_name(hn), "")) then
write (msg, '("Illegal surface index specified for surface ",i4," (index is ",i4,", must be less than ",i4,")")') hn-1, surface_index(hn), nsurfaces
call CCTK_PARAMWARN (msg)
goto 9999
end if
- sn = surface_index(hn) + 1
+ sn = sf_IdFromName(surface_index(hn), surface_name(hn)) + 1
! Import surface description
qlm_nghoststheta(hn) = nghoststheta(sn)
diff --git a/test/qlm-minkowski.par b/test/qlm-minkowski.par
index ddccd2c..ef9fb14 100644
--- a/test/qlm-minkowski.par
+++ b/test/qlm-minkowski.par
@@ -95,6 +95,7 @@ SphericalSurface::nsurfaces = 1
SphericalSurface::maxntheta = 39
SphericalSurface::maxnphi = 76
+SphericalSurface::name [0] = "Apparent Horizon"
SphericalSurface::ntheta [0] = 39
SphericalSurface::nphi [0] = 76
SphericalSurface::nghoststheta[0] = 2
@@ -155,7 +156,7 @@ QuasiLocalMeasures::verbose = yes
QuasiLocalMeasures::interpolator = "Lagrange polynomial interpolation"
QuasiLocalMeasures::interpolator_options = "order=4"
QuasiLocalMeasures::spatial_order = 4
-QuasiLocalMeasures::surface_index[0] = 0
+QuasiLocalMeasures::surface_name[0] = "Apparent Horizon"