1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
|
/*@@
@file AHFinder_mask.F
@date November 1998
@author Miguel Alcubierre
@desc
Set up horizon mask for black hole excision.
@enddesc
@@*/
#include "cctk.h"
#include "cctk_Parameters.h"
#include "cctk_Arguments.h"
subroutine AHFinder_mask(CCTK_FARGUMENTS,rhor)
! Set up the horizon mask. All this routine does is
! set up the value of the grid function "ahmask" to zero
! for points inside the horizon.
!
! It also finds out the radius of the largest sphere that
! fits inside the horizon and it centered in (xc,yc,zc).
! This is useful for simple excision.
use AHFinder_dat
implicit none
DECLARE_CCTK_FARGUMENTS
DECLARE_CCTK_PARAMETERS
CCTK_INT i,j,k
CCTK_REAL zero,one
CCTK_REAL aux,rhor
! *****************************
! *** DEFINE {zero,one} ***
! *****************************
zero = 0.0D0
one = 1.0D0
! *************************
! *** START ROUTINE ***
! *************************
! If inside of horizon (minus a bufferzone) set the mask to 0.
! To find out how big is the buffer zone inside theFor this
! I check the parameter ahf_maskshrink. This parameter
! says by which factor of the monopole term to shrink the
! mask with respect to the real horizon. Here I use the fact
! that the function ahfgrid is zero at the horizon, and is linear
! on the monopole term (which measures the overall scale).
!
! Notice that if I want a cubical mask, I first find out the radius
! of the largest sphere that fits inside the horizon and only after
! that I set the mask.
rhor = 1.0D10
call AHFinder_fun(CCTK_FARGUMENTS)
do k=1,nz
do j=1,ny
do i=1,nx
if (ahfgrid(i,j,k).ge.(ahf_maskshrink-one)*c0(0)) then
aux = dsqrt((x(i,j,k)-xc)**2 + (y(i,j,k)-yc)**2
. + (z(i,j,k)-xc)**2)
if (rhor.gt.aux) rhor = aux
else if (ahf_maskcube.eq.0) then
ahmask(i,j,k) = zero
end if
end do
end do
end do
! ************************
! *** CUBICAL MASK ***
! ************************
! If we want a cubical mask, set the mask based
! on the value of rhor.
if (ahf_maskcube.ne.0) then
do k=1,nz
do j=1,ny
do i=1,nx
if ((x(i,j,k).gt.(xc-0.5D0*rhor)).and.
. (y(i,j,k).gt.(yc-0.5D0*rhor)).and.
. (z(i,j,k).gt.(zc-0.5D0*rhor)).and.
. (x(i,j,k).lt.(xc+0.5D0*rhor)).and.
. (y(i,j,k).lt.(yc+0.5D0*rhor)).and.
. (z(i,j,k).lt.(zc+0.5D0*rhor))) then
ahmask(i,j,k) = zero
end if
end do
end do
end do
end if
! ****************************************************
! *** APPLY SYMMETRY BOUNDARIES AND SYNCHRONIZE ***
! ****************************************************
call CCTK_SyncGroup(cctkGH,"ahfinder::ahfmask")
call CartSymBCGroup(ierr,cctkGH,"ahfinder::ahfmask")
! ***************
! *** END ***
! ***************
end subroutine AHFinder_mask
|