aboutsummaryrefslogtreecommitdiff
path: root/src/GRHydro_Boundaries.F90
blob: 1c97c2c422914c71444ca026f51421260bc9b38b (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
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
 /*@@
   @file      GRHydro_Boundaries.F90
   @date      Sat Jan 26 01:01:14 2002
   @author    
   @desc 
   The two routines for dealing with boundary conditions.
   @enddesc 
 @@*/

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

#include "util_Table.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)
#define sx(i,j,k) scon(i,j,k,1)
#define sy(i,j,k) scon(i,j,k,2)
#define sz(i,j,k) scon(i,j,k,3)

 /*@@
   @routine    GRHydro_InitSymBound
   @date       Sat Jan 26 01:03:04 2002
   @author     Ian Hawke
   @desc 
   Sets up the symmetries at the boundaries of the hydrodynamical variables.
   @enddesc 
   @calls     
   @calledby   
   @history 
   Direct translation of routines from GR3D, GRAstro_Hydro, 
   written by Mark Miller, or WaveToy routines, or...
   @endhistory 

@@*/

subroutine GRHydro_InitSymBound(CCTK_ARGUMENTS)

  implicit none
  
  DECLARE_CCTK_ARGUMENTS
  DECLARE_CCTK_PARAMETERS
  DECLARE_CCTK_FUNCTIONS

  integer, dimension(3) :: sym
  integer :: ierr = 0
  integer :: itracer

  character(len=100) tracername
  character(len=100) tracerindex

  sym(1) = 1
  sym(2) = 1
  sym(3) = 1
  
  call SetCartSymVN(ierr, cctkGH, sym, "HydroBase::rho")
  call SetCartSymVN(ierr, cctkGH, sym, "HydroBase::press")
  call SetCartSymVN(ierr, cctkGH, sym, "GRHydro::dens")
  call SetCartSymVN(ierr, cctkGH, sym, "GRHydro::tau")
  call SetCartSymVN(ierr, cctkGH, sym, "GRHydro::w_lorentz")
  call SetCartSymVN(ierr, cctkGH, sym, "HydroBase::eps")
  call SetCartSymVN(ierr, cctkGH, sym, "GRHydro::GRHydro_C2P_failed")

!!$ handle multiple tracer variables
  if(evolve_tracer.ne.0) then
     call SetCartSymGN(ierr, cctkGH, sym, "GRHydro::GRHydro_tracers")
     call SetCartSymGN(ierr, cctkGH, sym, "GRHydro::GRHydro_cons_tracers")
  endif 

  if(evolve_y_e.ne.0) then
     call SetCartSymGN(ierr, cctkGH, sym, "HydroBase::Y_e")
     call SetCartSymGN(ierr, cctkGH, sym, "GRHydro::Y_e_con")
  endif

  if(evolve_temper.ne.0) then
     call SetCartSymGN(ierr, cctkGH, sym, "HydroBase::temperature")
     call SetCartSymGN(ierr, cctkGH, sym, "HydroBase::entropy")
  endif
  
  sym(1) = -1
  sym(2) = 1
  sym(3) = 1
  
  call SetCartSymVN(ierr, cctkGH, sym, "HydroBase::vel[0]")
  call SetCartSymVN(ierr, cctkGH, sym, "GRHydro::scon[0]")
  
  sym(1) = 1
  sym(2) = -1
  sym(3) = 1
  
  call SetCartSymVN(ierr, cctkGH, sym, "HydroBase::vel[1]")
  call SetCartSymVN(ierr, cctkGH, sym, "GRHydro::scon[1]")
  
  sym(1) = 1
  sym(2) = 1
  sym(3) = -1
  
  call SetCartSymVN(ierr, cctkGH, sym, "HydroBase::vel[2]")
  call SetCartSymVN(ierr, cctkGH, sym, "GRHydro::scon[2]")
  
end subroutine GRHydro_InitSymBound

 /*@@
   @routine    GRHydro_Boundaries
   @date       Sat Jan 26 01:04:04 2002
   @author     
   @desc 
   Calls the appropriate boundary routines
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/

subroutine GRHydro_Boundaries(CCTK_ARGUMENTS)
  
  implicit none
  
  DECLARE_CCTK_ARGUMENTS
  DECLARE_CCTK_PARAMETERS
  DECLARE_CCTK_FUNCTIONS
  
  integer, dimension(3) :: sw
  integer :: ierr
  CCTK_REAL :: pi = 3.141569d0
  integer :: i,j,k

  CCTK_INT, parameter :: faces=CCTK_ALL_FACES
  CCTK_INT, parameter :: ione=1
  
  sw = GRHydro_stencil

!!$Flat boundaries if required  

  if (CCTK_EQUALS(bound,"flat")) then
    ierr = Boundary_SelectGroupForBC(cctkGH, faces, GRHydro_stencil, -ione, &
         "GRHydro::dens", "Flat")
    ierr = Boundary_SelectGroupForBC(cctkGH, faces, GRHydro_stencil, -ione, &
         "GRHydro::tau", "Flat")
    ierr = Boundary_SelectGroupForBC(cctkGH, faces, GRHydro_stencil, -ione, &
         "GRHydro::scon", "Flat")
    ierr = Boundary_SelectGroupForBC(cctkGH, faces, GRHydro_stencil, -ione, &
         "GRHydro::w_lorentz", "Flat")
    ierr = Boundary_SelectGroupForBC(cctkGH, faces, GRHydro_stencil, -ione, &
         "HydroBase::rho", "Flat")
    ierr = Boundary_SelectGroupForBC(cctkGH, faces, GRHydro_stencil, -ione, &
         "HydroBase::press", "Flat")
    ierr = Boundary_SelectGroupForBC(cctkGH, faces, GRHydro_stencil, -ione, &
         "HydroBase::eps", "Flat")
    ierr = Boundary_SelectGroupForBC(cctkGH, faces, GRHydro_stencil, -ione, &
         "HydroBase::vel", "Flat")

    if(evolve_tracer.ne.0) then 
      ierr = Boundary_SelectGroupForBC(cctkGH, faces, GRHydro_stencil, -ione, &
            "GRHydro::GRHydro_tracers", "Flat")
      ierr = Boundary_SelectGroupForBC(cctkGH, faces, GRHydro_stencil, -ione, &
            "GRHydro::GRHydro_cons_tracers", "Flat")
    endif

    if(evolve_y_e.ne.0) then
      ierr = Boundary_SelectGroupForBC(cctkGH, faces, GRHydro_stencil, -ione, &
            "HydroBase::Y_e", "Flat")
      ierr = Boundary_SelectGroupForBC(cctkGH, faces, GRHydro_stencil, -ione, &
            "GRHydro::Y_e_con", "Flat")
    endif

    if(evolve_temper.ne.0) then
      ierr = Boundary_SelectGroupForBC(cctkGH, faces, GRHydro_stencil, -ione, &
            "HydroBase::temperature", "Flat")
      ierr = Boundary_SelectGroupForBC(cctkGH, faces, GRHydro_stencil, -ione, &
            "HydroBase::entropy", "Flat")
    endif


  endif

  if (CCTK_EQUALS(bound,"none")) then
    ierr = Boundary_SelectGroupForBC(cctkGH, faces, GRHydro_stencil, -ione, &
         "GRHydro::dens", "None")
    ierr = Boundary_SelectGroupForBC(cctkGH, faces, GRHydro_stencil, -ione, &
         "GRHydro::tau", "None")
    ierr = Boundary_SelectGroupForBC(cctkGH, faces, GRHydro_stencil, -ione, &
         "GRHydro::scon", "None")
    ierr = Boundary_SelectGroupForBC(cctkGH, faces, GRHydro_stencil, -ione, &
         "GRHydro::w_lorentz", "None")
    ierr = Boundary_SelectGroupForBC(cctkGH, faces, GRHydro_stencil, -ione, &
         "HydroBase::rho", "None")
    ierr = Boundary_SelectGroupForBC(cctkGH, faces, GRHydro_stencil, -ione, &
         "HydroBase::press", "None")
    ierr = Boundary_SelectGroupForBC(cctkGH, faces, GRHydro_stencil, -ione, &
         "HydroBase::eps", "None")
    ierr = Boundary_SelectGroupForBC(cctkGH, faces, GRHydro_stencil, -ione, &
         "HydroBase::vel", "None")

    if(evolve_tracer.ne.0) then 
      ierr = Boundary_SelectGroupForBC(cctkGH, faces, GRHydro_stencil, -ione, &
            "GRHydro::GRHydro_tracers", "None")
      ierr = Boundary_SelectGroupForBC(cctkGH, faces, GRHydro_stencil, -ione, &
            "GRHydro::GRHydro_cons_tracers", "None")
    endif

    if(evolve_y_e.ne.0) then
      ierr = Boundary_SelectGroupForBC(cctkGH, faces, GRHydro_stencil, -ione, &
            "HydroBase::Y_e", "None")
      ierr = Boundary_SelectGroupForBC(cctkGH, faces, GRHydro_stencil, -ione, &
            "GRHydro::Y_e_con", "None")
    endif

    if(evolve_temper.ne.0) then
      ierr = Boundary_SelectGroupForBC(cctkGH, faces, GRHydro_stencil, -ione, &
            "HydroBase::temperature", "None")
      ierr = Boundary_SelectGroupForBC(cctkGH, faces, GRHydro_stencil, -ione, &
            "HydroBase::entropy", "None")
    endif

  end if

  if (CCTK_EQUALS(bound,"scalar")) then
    call CCTK_WARN(0, "Until somebody uses this I see no reason to support it")
  end if

  if (ierr < 0) call CCTK_WARN(0, "problems with applying the chosen boundary condition")

end subroutine GRHydro_Boundaries


 /*@@
   @routine    GRHydro_OutflowBoundaries
   @date       Tue May 17 16:58:06 2005
   @author     Ian Hawke
   @desc 
   Set outflow boundaries over only part of the domain.
   This is designed to be used with GZPatchSystem.
   @enddesc 
   @calls     
   @calledby   
   @history 
 
   @endhistory 

@@*/

subroutine GRHydro_OutflowBoundaries(CCTK_ARGUMENTS)
  
  implicit none
  
  DECLARE_CCTK_ARGUMENTS
  DECLARE_CCTK_PARAMETERS
  DECLARE_CCTK_FUNCTIONS
  
  integer, dimension(3) :: sw
  integer :: ierr
  integer :: i,j,k
  
  CCTK_REAL, dimension(3) :: posn

  CCTK_REAL :: det,psi4pt

  sw = GRHydro_stencil

  if (r(1,1,1) < r(1,1,cctk_lsh(3))) then

    if (cctk_bbox(6) .ne. 0) then

      do k = cctk_lsh(3) - sw(3), cctk_lsh(3)
        do j = 1, cctk_lsh(2)
          do i = 1, cctk_lsh(1)
            
            posn(1) = x(i,j,k)
            posn(2) = y(i,j,k)
            posn(3) = z(i,j,k)
            
            if (dot_product(outflowboundary_normal, posn) > 0.d0) then
              
              rho(i,j,k)        = rho(i,j,k-1)
              velx(i,j,k)       = velx(i,j,k-1)
              vely(i,j,k)       = vely(i,j,k-1)
              velz(i,j,k)       = velz(i,j,k-1)
              eps(i,j,k)        = eps(i,j,k-1)
              press(i,j,k)      = press(i,j,k-1)
              w_lorentz(i,j,k)  = w_lorentz(i,j,k-1)
              
              psi4pt = 1.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 prim2con(GRHydro_eos_handle,psi4pt*gxx(i,j,k),&
                   psi4pt*gxy(i,j,k),psi4pt*gxz(i,j,k),&
                   psi4pt*gyy(i,j,k),psi4pt*gyz(i,j,k),psi4pt*gzz(i,j,k),&
                   det, dens(i,j,k),sx(i,j,k),sy(i,j,k),sz(i,j,k),&
                   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))
              
            end if
          
          end do
        end do
      end do
      
    end if
  
  else

    if (cctk_bbox(5) .ne. 0) then

      do k = sw(3), 1, -1
        do j = 1, cctk_lsh(2)
          do i = 1, cctk_lsh(1)
            
            posn(1) = x(i,j,k)
            posn(2) = y(i,j,k)
            posn(3) = z(i,j,k)
            
            if (dot_product(outflowboundary_normal, posn) > 0.d0) then
              
              rho(i,j,k)        = rho(i,j,k+1)
              velx(i,j,k)       = velx(i,j,k+1)
              vely(i,j,k)       = vely(i,j,k+1)
              velz(i,j,k)       = velz(i,j,k+1)
              eps(i,j,k)        = eps(i,j,k+1)
              press(i,j,k)      = press(i,j,k+1)
              w_lorentz(i,j,k)  = w_lorentz(i,j,k+1)
              
              psi4pt = 1.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 prim2con(GRHydro_eos_handle,psi4pt*gxx(i,j,k),&
                   psi4pt*gxy(i,j,k),psi4pt*gxz(i,j,k),&
                   psi4pt*gyy(i,j,k),psi4pt*gyz(i,j,k),psi4pt*gzz(i,j,k),&
                   det, dens(i,j,k),sx(i,j,k),sy(i,j,k),sz(i,j,k),&
                   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))
              
            end if
          
          end do
        end do
      end do
      
    end if

  end if
  
end subroutine GRHydro_OutflowBoundaries