From 6687c44ec86e4aaf48b14a5ee9a4d2638b515d42 Mon Sep 17 00:00:00 2001 From: diener Date: Fri, 1 Aug 2003 14:34:43 +0000 Subject: Now uses SetDriftCorrectPosition to set the Driftcorrect centroid. git-svn-id: http://svn.einsteintoolkit.org/cactus/EinsteinAnalysis/AHFinder/trunk@352 89daf98e-ef62-4674-b946-b8ff9de2216c --- src/AHFinder_gau.F | 68 +++++++++++++++++++++++++++++++++++++++++------------- 1 file changed, 52 insertions(+), 16 deletions(-) diff --git a/src/AHFinder_gau.F b/src/AHFinder_gau.F index 0e97711..5de6103 100644 --- a/src/AHFinder_gau.F +++ b/src/AHFinder_gau.F @@ -68,6 +68,8 @@ CCTK_REAL, allocatable, dimension(:,:) :: txx,tyy,tzz,txy,txz,tyz CCTK_REAL, allocatable, dimension(:,:) :: g11,g22,g12 + CCTK_REAL, dimension(3) :: ahf_centroid, rahf_centroid + CCTK_REAL :: m1p, m2p, x1p, x2p, y1p, y2p, z1p, z2p CCTK_REAL :: r12, r22, sigma CCTK_REAL, dimension(nx,ny,nz) :: fac @@ -76,6 +78,7 @@ logical :: str_comp logical :: gaussf_exists external str_comp +! CCTK_INT :: CCTK_IsFunctionAliased character(len=32) :: tmpstr @@ -605,24 +608,15 @@ avgy = avgy/dble(npoints) avgz = avgz/dble(npoints) - if ((.not.find3).or. - . (find3.and.(mfind-1.eq.horizon_centroid))) then - - ahf_centroid_x = avgx - ahf_centroid_y = avgy - -! Take symmetries into account. + ahf_centroid(1) = avgx + ahf_centroid(2) = avgy + ahf_centroid(3) = avgz - if (refx) ahf_centroid_x = xc - if (refy) ahf_centroid_y = yc - - if (verbose) then - write(*,*) - write(6,*) 'AHFinder_gau: ahf_centroid_x = ', ahf_centroid_x - write(6,*) 'AHFinder_gau: ahf_centroid_y = ', ahf_centroid_y - end if +! Take symmetries into account. - end if + if (refx) ahf_centroid(1) = xc + if (refy) ahf_centroid(2) = yc + if (refz) ahf_centroid(3) = zc ! Other processors. @@ -632,6 +626,9 @@ ya = half*(ymx+ymn) za = half*(zmx+zmn) + ahf_centroid(1) = zero + ahf_centroid(2) = zero + ahf_centroid(3) = zero end if @@ -655,6 +652,45 @@ end if error2 = rerror + call CCTK_ReduceLocArrayToArray1D(ierror,cctkGH,-1,sum_handle, + . ahf_centroid,rahf_centroid,3,CCTK_VARIABLE_REAL) + + if (ierror.ne.0) then + call CCTK_WARN(1,"Reduction failed!") + end if + + if ( ( which_horizon_to_announce_centroid .gt. 0 ) .and. + . ( which_horizon_to_announce_centroid .eq. mfind ) ) then + + call CCTK_IsFunctionAliased(ierror, + . "SetDriftCorrectPosition") + if ( ierror .eq. 1 ) then + call CCTK_INFO("Announcing centroid to DriftCorrect") + call SetDriftCorrectPosition ( cctkGH, ahf_centroid(1), + . ahf_centroid(2), + . ahf_centroid(3) ) + end if + end if + + if ( ( which_horizon_to_output_centroid .gt. 0 ) .and. + . ( which_horizon_to_output_centroid .eq. mfind ) ) then + + ahf_centroid_x = rahf_centroid(1) + ahf_centroid_y = rahf_centroid(2) + ahf_centroid_z = rahf_centroid(3) + + if (verbose) then + write(*,*) + write(6,*) 'AHFinder_gau: ahf_centroid_x = ', + . ahf_centroid_x + write(6,*) 'AHFinder_gau: ahf_centroid_y = ', + . ahf_centroid_y + write(6,*) 'AHFinder_gau: ahf_centroid_z = ', + . ahf_centroid_z + end if + end if + + ! If there was an error then return from the subroutine ! (but remember to deallocate the arrays first!). -- cgit v1.2.3