aboutsummaryrefslogtreecommitdiff
path: root/src/EHFinder_Init.F90
diff options
context:
space:
mode:
authordiener <diener@2a26948c-0e4f-0410-aee8-f1d3e353619c>2002-07-16 11:31:41 +0000
committerdiener <diener@2a26948c-0e4f-0410-aee8-f1d3e353619c>2002-07-16 11:31:41 +0000
commit1789cb89e8d4354346b557c47b9e35d397d0db2b (patch)
tree5f5008766cd55a6da9660ed2b685326e3a06bfd5 /src/EHFinder_Init.F90
parente1091474db16b640eb236ba021ec10bf3acc73e0 (diff)
Defines ex_value to be the value of f in excised points to avoid complications
when making contour plots or isosurfaces. Modifies cctk_time to initially being equal to the final time in the numerical evolution. git-svn-id: http://svn.einsteintoolkit.org/cactus/EinsteinAnalysis/EHFinder/trunk@33 2a26948c-0e4f-0410-aee8-f1d3e353619c
Diffstat (limited to 'src/EHFinder_Init.F90')
-rw-r--r--src/EHFinder_Init.F908
1 files changed, 8 insertions, 0 deletions
diff --git a/src/EHFinder_Init.F90 b/src/EHFinder_Init.F90
index 4535268..2e4ebad 100644
--- a/src/EHFinder_Init.F90
+++ b/src/EHFinder_Init.F90
@@ -17,8 +17,16 @@ subroutine EHFinder_Init(CCTK_ARGUMENTS)
! CCTK_INT, dimension(0:63) :: eh_check
CCTK_INT :: i, j, k
+ CCTK_REAL :: last_time
+
+ ex_value = - ( one + shell_width ) * maxval(cctk_delta_space)
current_iteration = last_iteration_number
+
+ last_time = abs(cctk_delta_time) * last_iteration_number
+
+ cctk_time = last_time
+
if ( CCTK_EQUALS( initial_f, 'sphere' ) ) then
f = sqrt(x**2+y**2+z**2) - initial_rad ! + &
! 0.5d0 * sin ( 5d0 * ( sqrt(x**2+y**2+z**2) - initial_rad ) )