aboutsummaryrefslogtreecommitdiff
path: root/src/GRHydro_UpdateMask.F90
blob: d7cdff0d7b3ddd9edca493d4c06241147ead0f17 (plain)
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
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
 /*@@
   @file      GRHydro_UpdateMask.F90
   @date      Wed Mar 13 14:18:38 2002
   @author    
   @desc 
   Alter the update terms if inside the atmosphere or excision region
   @enddesc 
 @@*/

#include "cctk.h"
#include "cctk_Parameters.h"
#include "cctk_Arguments.h"

#include "GRHydro_Macros.h"
#include "SpaceMask.h"

#define velx(i,j,k) vel(i,j,k,1)
#define vely(i,j,k) vel(i,j,k,2)
#define velz(i,j,k) vel(i,j,k,3)

subroutine GRHydroUpdateAtmosphereMask(CCTK_ARGUMENTS)

  implicit none
  
  DECLARE_CCTK_ARGUMENTS
  DECLARE_CCTK_PARAMETERS

  CCTK_INT :: i,j,k
  CCTK_REAL :: frac

  CCTK_INT :: type_bits, atmosphere

  call SpaceMask_GetTypeBits(type_bits, "Hydro_Atmosphere")
  call SpaceMask_GetStateBits(atmosphere, "Hydro_Atmosphere", "in_atmosphere")

  frac = CCTK_DELTA_TIME

  !$OMP PARALLEL DO PRIVATE(j,i)
  do k = 1, cctk_lsh(3)
    do j = 1, cctk_lsh(2)
      do i = 1, cctk_lsh(1)
        if ( GRHydro_enable_internal_excision /= 0 .and. (hydro_excision_mask(i,j,k) .ne. 0) .or. &
             (SpaceMask_CheckStateBitsF90(space_mask, i, j, k, \
                                          type_bits, atmosphere)) .or. &
             (tau(i,j,k) + frac * taurhs(i,j,k) .le. 0.d0) .or. &
             (dens(i,j,k) + frac * densrhs(i,j,k) .le. 0.d0) ) then
          densrhs(i,j,k) = 0.0d0
          srhs(i,j,k,:)   = 0.0d0
          taurhs(i,j,k)  = 0.0d0
          atmosphere_mask(i,j,k) = 1

          SpaceMask_SetStateBitsF90(space_mask, i, j, k, type_bits, atmosphere)
        end if
      end do
    end do
  end do
  !$OMP END PARALLEL DO

end subroutine GRHydroUpdateAtmosphereMask

 /*@@
   @routine    GRHydro_SetupMask
   @date       Thu Jun 20 13:27:28 2002
   @author     Ian Hawke
   @desc 
   Initialize the mask to be zero.
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/

subroutine GRHydro_SetupMask(CCTK_ARGUMENTS)

  implicit none

  DECLARE_CCTK_ARGUMENTS
  DECLARE_CCTK_PARAMETERS

  CCTK_INT :: type_bits, not_atmosphere, i, j, k

  ! Initialize all rhs variables and the mask.
  ! The former vars need to be initialized since there is
  ! no rhs computation in CCTK_INITIAL or POSTINITIAL.

  densrhs = 0.0d0
  taurhs  = 0.0d0 
  srhs    = 0.0d0

  if (evolve_tracer .ne. 0) then
     cons_tracerrhs = 0.0d0
  endif

  atmosphere_mask = 0
  
  call SpaceMask_GetTypeBits(type_bits, "Hydro_Atmosphere")
  call SpaceMask_GetStateBits(not_atmosphere, &
       &"Hydro_Atmosphere", "not_in_atmosphere")

  do k = 1, cctk_lsh(3)
    do j = 1, cctk_lsh(2)
      do i = 1, cctk_lsh(1)
  
        SpaceMask_SetStateBitsF90(space_mask, i, j, k, type_bits, not_atmosphere)

      end do
    end do
  end do

  call CCTK_INFO("Setting up the atmosphere mask: all points are not_atmosphere")
  
end subroutine GRHydro_SetupMask

 /*@@
   @routine    GRHydro_InitAtmosMask
   @date       Thu Jun 20 13:27:28 2002
   @author     Ian Hawke
   @desc 
   Initialize the mask based on rho_min. This is used only if wk_atmosphere=yes.
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/

subroutine GRHydro_InitAtmosMask(CCTK_ARGUMENTS)

  implicit none

  DECLARE_CCTK_ARGUMENTS

  CCTK_INT :: type_bits, atmosphere, i, j, k

  call SpaceMask_GetTypeBits(type_bits, "Hydro_Atmosphere")
  call SpaceMask_GetStateBits(atmosphere, &
       &"Hydro_Atmosphere", "in_atmosphere")

  do k = 1, cctk_lsh(3)
    do j = 1, cctk_lsh(2)
      do i = 1, cctk_lsh(1)

        if (rho(i,j,k) .le. GRHydro_rho_min) then
  
          SpaceMask_SetStateBitsF90(space_mask, i, j, k, type_bits, atmosphere)
          atmosphere_mask(i,j,k) = 1

        end if
        
      end do
    end do
  end do

  call CCTK_INFO("Setting up the atmosphere mask: points with rho<rho_min are set to be atmosphere")
  
end subroutine GRHydro_InitAtmosMask

 /*@@
   @routine    GRHydro_AtmosphereReset
   @date       Thu Jun 20 13:30:51 2002
   @author     Ian Hawke
   @desc 
   After MoL has evolved, if a point is supposed to be reset then do so.
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/

subroutine GRHydro_AtmosphereReset(CCTK_ARGUMENTS)

  implicit none

  DECLARE_CCTK_ARGUMENTS
  DECLARE_CCTK_PARAMETERS

  CCTK_INT :: i, j, k
  CCTK_REAL :: det, psi4pt

  CCTK_INT :: type_bits, atmosphere, not_atmosphere
  
  call SpaceMask_GetTypeBits(type_bits, "Hydro_Atmosphere")
  call SpaceMask_GetStateBits(atmosphere, "Hydro_Atmosphere",&
                              "in_atmosphere")
  call SpaceMask_GetStateBits(not_atmosphere, "Hydro_Atmosphere",&
                              "not_in_atmosphere")

  do k = 1, cctk_lsh(3)
    do j = 1, cctk_lsh(2)
      do i = 1, cctk_lsh(1)
        
        if ( (atmosphere_mask(i, j, k) .eq. 1) &
             &.or. (SpaceMask_CheckStateBitsF90(space_mask,i, j, k, type_bits,\
                    atmosphere)) &
             &) then

!!$          write(*,*) 'Resetting at ',i,j,k, atmosphere_mask(i, j, k), &
!!$             & (SpaceMask_CheckStateBitsF90(space_mask,i, j, k, type_bits,\
!!$                    atmosphere)) 

          rho(i,j,k) = GRHydro_rho_min
          velx(i,j,k) = 0.0d0
          vely(i,j,k) = 0.0d0
          velz(i,j,k) = 0.0d0
          det = SPATIAL_DETERMINANT(gxx(i,j,k), gxy(i,j,k), gxz(i,j,k), \
               gyy(i,j,k), gyz(i,j,k), gzz(i,j,k))
          call prim2conpolytype(GRHydro_polytrope_handle, &
               gxx(i,j,k), gxy(i,j,k), gxz(i,j,k), &
               gyy(i,j,k), gyz(i,j,k), gzz(i,j,k), det, &
               dens(i,j,k), scon(i,j,k,1), scon(i,j,k,2), scon(i,j,k,3), &
               tau(i,j,k), rho(i,j,k), velx(i,j,k), vely(i,j,k), &
               velz(i,j,k), eps(i,j,k), press(i,j,k), w_lorentz(i,j,k))
          if (wk_atmosphere .eq. 0) then
            atmosphere_mask(i, j, k) = 0
            SpaceMask_SetStateBitsF90(space_mask, i, j, k, type_bits,\
                                      not_atmosphere)
          end if

        end if

      end do
    end do
  end do

!!$  call GRHydro_Boundaries(CCTK_PASS_FTOF)

end subroutine GRHydro_AtmosphereReset