-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathbc.F90
384 lines (325 loc) · 12.8 KB
/
bc.F90
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
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
module bc
#include "definition.h"
use DMR
use grid_data
use sim_data, only : sim_xBC, sim_yBC, sim_bcTypex, sim_bcTypey, sim_bcTypez, sim_reconmultiD
implicit none
contains
subroutine bc_init(bc_char, bc_int)
implicit none
character(len=MAX_STRING_LENGTH), intent(IN ) :: bc_char
integer, intent(OUT) :: bc_int
if (bc_char == 'periodic') then
bc_int = PERIODIC
elseif (bc_char == 'outflow') then
bc_int = OUTFLOW
elseif (bc_char == 'DMR') then
bc_int = DBLMCH
elseif (bc_char == 'inflow') then
bc_int = INFLOW
elseif (bc_char == 'reflect') then
bc_int = REFLECT
end if
end subroutine bc_init
!!$ subroutine bc_apply(V)
!!$ implicit none
!!$ real, dimension(NUMB_VAR,gr_imax(XDIM),gr_imax(YDIM)), intent(INOUT) :: V
!!$ return
!!$ end subroutine bc_apply
subroutine bc_apply(V)
use block_data
implicit none
real, dimension(NUMB_VAR,gr_imax(XDIM),gr_imax(YDIM),gr_imax(ZDIM)), intent(INOUT) :: V
integer :: dest, gc, i,j,k
real, dimension(NUMB_VAR,gr_ngc, gr_ny ,gr_nz ) :: loc_buffL, loc_buffR
real, dimension(NUMB_VAR,gr_nx , gr_ngc,gr_nz ) :: loc_buffT, loc_buffB
real, dimension(NUMB_VAR,gr_nx , gr_ny ,gr_ngc) :: loc_buffU, loc_buffD
!BCS: 1: left x L
! 2: right y T
! 3: right x R
! 4: left y B
! 5: left z D
! 6: right z U
!!$ !load halo data into image buffers
bl_buffL(:, :, :, :) = &
V(:, gr_ibeg(XDIM):gr_ibeg(XDIM)+gr_ngc-1, gr_ibeg(YDIM):gr_iend(YDIM),gr_ibeg(ZDIM):gr_iend(ZDIM))
bl_buffR(:, :, :, :) = &
V(:, gr_iend(XDIM)-gr_ngc+1:gr_iend(XDIM), gr_ibeg(YDIM):gr_iend(YDIM),gr_ibeg(ZDIM):gr_iend(ZDIM))
bl_buffB(:, :, :, :) = &
V(:, gr_ibeg(XDIM):gr_iend(XDIM), gr_ibeg(YDIM):gr_ibeg(YDIM)+gr_ngc-1,gr_ibeg(ZDIM):gr_iend(ZDIM))
bl_buffT(:, :, :, :) = &
V(:, gr_ibeg(XDIM):gr_iend(XDIM), gr_iend(YDIM)-gr_ngc+1:gr_iend(YDIM),gr_ibeg(ZDIM):gr_iend(ZDIM))
bl_buffD(:, :, :, :) = &
V(:, gr_ibeg(XDIM):gr_iend(XDIM), gr_ibeg(YDIM):gr_iend(YDIM),gr_ibeg(ZDIM):gr_ibeg(ZDIM)+gr_ngc-1)
bl_buffU(:, :, :, :) = &
V(:, gr_ibeg(XDIM):gr_iend(XDIM), gr_ibeg(YDIM):gr_iend(YDIM),gr_iend(ZDIM)-gr_ngc+1:gr_iend(ZDIM))
Sync All
!do halo exchanges along non-domain boundaries
dest = bl_BC(1)
if (dest > 0) then
loc_buffL(:,:,:,:) = bl_buffR(:,:,:,:)[dest]
end if
dest = bl_BC(2)
if (dest > 0) then
loc_buffT(:,:,:,:) = bl_buffB(:,:,:,:)[dest]
end if
dest = bl_BC(3)
if (dest > 0) then
loc_buffR(:,:,:,:) = bl_buffL(:,:,:,:)[dest]
end if
dest = bl_BC(4)
if (dest > 0) then
loc_buffB(:,:,:,:) = bl_buffT(:,:,:,:)[dest]
end if
dest = bl_BC(5)
if (dest > 0) then
loc_buffD(:,:,:,:) = bl_buffU(:,:,:,:)[dest]
end if
dest = bl_BC(6)
if (dest > 0) then
loc_buffU(:,:,:,:) = bl_buffD(:,:,:,:)[dest]
end if
!now we apply domain boundary conditions
if (sim_bcTypex == 'periodic') then
!for periodic we can just handle exchange here
if (bl_i == 1 ) then
dest = bl_grid(bl_iProcs, bl_j, bl_k)
loc_buffL(:,:,:,:) = bl_buffR(:,:,:,:)[dest]
end if
if (bl_i == bl_iProcs) then
dest = bl_grid(1 , bl_j, bl_k)
loc_buffR(:,:,:,:) = bl_buffL(:,:,:,:)[dest]
end if
elseif (sim_bcTypex == 'reflect') then
if (bl_i == 1) then
loc_buffL(:,:,:,:) = V(:,gr_ibeg(XDIM)+gr_ngc:gr_ibeg(XDIM):-1, &
gr_ibeg(YDIM):gr_iend(YDIM),gr_ibeg(ZDIM):gr_iend(ZDIM))
loc_buffL(VELX_VAR,:,:,:) = -loc_buffL(VELX_VAR,:,:,:)
end if
if (bl_i == bl_iProcs) then
loc_buffR(:,:,:,:) = V(:,gr_iend(XDIM)-gr_ngc:gr_iend(XDIM):-1, &
gr_ibeg(YDIM):gr_iend(YDIM),gr_ibeg(ZDIM):gr_iend(ZDIM))
loc_buffR(VELX_VAR,:,:,:) = -loc_buffR(VELX_VAR,:,:,:)
end if
elseif (sim_bcTypex == 'inflow') then
if (bl_i == 1) then
!let's leave as initial conditions on the inflow side
loc_buffL(:,:,:,:) = V(:,gr_i0 (XDIM):gr_ibeg(XDIM)-1, &
gr_ibeg(YDIM):gr_iend(YDIM) , gr_ibeg(ZDIM):gr_iend(ZDIM))
end if
if (bl_i == bl_iProcs) then
!outflow on the other side
do i = gr_iend(XDIM)+1, gr_imax(XDIM)
loc_buffR(:,i-gr_iend(XDIM),:,:) = V(:,gr_iend(XDIM), &
gr_ibeg(YDIM):gr_iend(YDIM) , gr_ibeg(ZDIM):gr_iend(ZDIM))
end do
end if
elseif (sim_bcTypex == 'outflow') then
if (bl_i == 1 ) then
do i = gr_i0(XDIM), gr_ibeg(XDIM)-1
loc_buffL(:,i,:,:) = V(:,gr_ibeg(XDIM), &
gr_ibeg(YDIM):gr_iend(YDIM) , gr_ibeg(ZDIM):gr_iend(ZDIM))
end do
end if
if (bl_i == bl_iProcs) then
do i = gr_iend(XDIM)+1, gr_imax(XDIM)
loc_buffR(:,i-gr_iend(XDIM),:,:) = V(:,gr_iend(XDIM), &
gr_ibeg(YDIM):gr_iend(YDIM) , gr_ibeg(ZDIM):gr_iend(ZDIM))
end do
end if
end if
if (sim_bcTypey == 'periodic') then
!again for periodic we can handle exchange here
if (bl_j == 1) then
dest = bl_grid(bl_i, bl_jProcs, bl_k)
loc_buffB(:,:,:,:) = bl_buffT(:,:,:,:)[dest]
end if
if (bl_j == bl_jProcs) then
dest = bl_grid(bl_i, 1, bl_k)
loc_buffT(:,:,:,:) = bl_buffB(:,:,:,:)[dest]
end if
elseif (sim_bcTypey == 'outflow') then
if (bl_j == 1) then
do j = gr_i0(YDIM), gr_ibeg(YDIM)-1
loc_buffB(:,:,j,:) = V(:, gr_ibeg(XDIM):gr_iend(XDIM), gr_ibeg(YDIM), &
gr_ibeg(ZDIM):gr_iend(ZDIM) )
end do
end if
if (bl_j == bl_jProcs) then
do j = gr_iend(YDIM)+1, gr_imax(YDIM)
loc_buffT(:,:,j-gr_iend(YDIM),:) = V(:, gr_ibeg(XDIM):gr_iend(XDIM), gr_iend(YDIM), &
gr_ibeg(ZDIM):gr_iend(ZDIM) )
end do
end if
elseif(sim_bcTypey == 'reflect') then
if (bl_j == 1) then
do i = gr_ibeg(XDIM), gr_iend(XDIM)
loc_buffB(:,i-gr_ngc,:,:) = V(:,i,gr_ibeg(YDIM)+gr_ngc:gr_ibeg(YDIM):-1,gr_ibeg(ZDIM):gr_iend(ZDIM))
loc_buffB(VELY_VAR,i-gr_ngc,:,:) = -loc_buffB(VELY_VAR,i-gr_ngc,:,:)
end do
end if
if (bl_j == bl_jProcs) then
loc_buffT(:,:,:,:) = V(:,gr_ibeg(XDIM):gr_iend(XDIM), &
gr_iend(YDIM):gr_iend(YDIM)-gr_ngc:-1, gr_ibeg(ZDIM):gr_iend(ZDIM))
loc_buffT(VELY_VAR,:,:,:) = -loc_buffT(VELY_VAR,:,:,:)
end if
end if
if (sim_bcTypez == 'periodic') then
if (bl_k == 1) then
dest = bl_grid(bl_i, bl_j, bl_kProcs)
loc_buffD(:,:,:,:) = bl_buffU(:,:,:,:)[dest]
end if
if (bl_k == bl_kProcs) then
dest = bl_grid(bl_i, bl_j, 1)
loc_buffU(:,:,:,:) = bl_buffD(:,:,:,:)[dest]
end if
elseif (sim_bcTypez == 'outflow') then
if (bl_k == 1) then
do k = gr_i0(ZDIM), gr_ibeg(ZDIM)-1
loc_buffD(:,:,:,k) = V(:,gr_ibeg(XDIM):gr_iend(XDIM), &
gr_ibeg(YDIM):gr_iend(YDIM), &
gr_ibeg(ZDIM))
end do
end if
if (bl_k == bl_kProcs) then
do k = gr_iend(ZDIM)+1, gr_imax(ZDIM)
loc_buffU(:,:,:,k-gr_iend(ZDIM)) = V(:,gr_ibeg(XDIM):gr_iend(XDIM), &
gr_ibeg(YDIM):gr_iend(YDIM), &
gr_iend(ZDIM))
end do
end if
elseif (sim_bcTypez == 'reflect') then
if (bl_k == 1) then
loc_buffD(:,:,:,:) = V(:,gr_ibeg(XDIM):gr_iend(XDIM), gr_ibeg(YDIM):gr_iend(YDIM), &
gr_ibeg(ZDIM)+gr_ngc:gr_ibeg(ZDIM):-1)
loc_buffD(VELZ_VAR,:,:,:) = - loc_buffD(VELZ_VAR,:,:,:)
end if
if (bl_k == bl_kProcs) then
loc_buffU(:,:,:,:) = V(:,gr_ibeg(XDIM):gr_iend(XDIM), gr_ibeg(YDIM):gr_iend(YDIM), &
gr_iend(ZDIM)-gr_ngc:gr_iend(ZDIM):-1)
loc_buffU(VELZ_VAR,:,:,:) = - loc_buffU(VELZ_VAR,:,:,:)
end if
end if
!fill halo data
!do this guard cell by guard cell so that it fills int he same order
!x GC's
V(:, gr_i0(XDIM):gr_ibeg(XDIM)-1, gr_ibeg(YDIM):gr_iend(YDIM),gr_ibeg(ZDIM):gr_iend(ZDIM)) = &
loc_buffL(:,:,:,:)
V(:, gr_iend(XDIM)+1:gr_imax(XDIM), gr_ibeg(YDIM):gr_iend(YDIM),gr_ibeg(ZDIM):gr_iend(ZDIM)) = &
loc_buffR(:,:,:,:)
!y GC's
V(:, gr_ibeg(XDIM):gr_iend(XDIM), gr_iend(YDIM)+1:gr_imax(YDIM),gr_ibeg(ZDIM):gr_iend(ZDIM)) = &
loc_buffT(:,:,:,:)
V(:, gr_ibeg(XDIM):gr_iend(XDIM), gr_i0(YDIM):gr_ibeg(YDIM)-1,gr_ibeg(ZDIM):gr_iend(ZDIM)) = &
loc_buffB(:,:,:,:)
!z GC's
V(:,gr_ibeg(XDIM):gr_iend(XDIM), gr_ibeg(YDIM):gr_iend(YDIM), gr_i0(ZDIM):gr_ibeg(ZDIM)-1) = &
loc_buffD(:,:,:,:)
V(:,gr_ibeg(XDIM):gr_iend(XDIM), gr_ibeg(YDIM):gr_iend(YDIM), gr_iend(ZDIM)+1:gr_imax(ZDIM)) = &
loc_buffU(:,:,:,:)
! if (sim_reconMultiD) call bc_corners(V)
end subroutine bc_apply
subroutine bc_outflowX(V, loc_buffL, loc_buffR)
implicit none
real, dimension(NUMB_VAR,gr_imax(XDIM),gr_imax(YDIM)), intent(INOUT) :: V
real, dimension(NUMB_VAR, gr_ngc, gr_ny), intent(INOUT) :: loc_buffL, loc_buffR
real, dimension(NUMB_VAR) :: Vbeg, Vend
integer :: i,j
!x-region
!!$ loc_buffL(:,:,:) = V(:, gr_iend(XDIM):gr_imax(XDIM)-1)
!!$ do i = 1, gr_ngc
!!$ loc_buffL
!!$ end do
!!$ do j = gr_ibeg(YDIM),gr_iend(YDIM)!1,gr_ny
!!$ do i = 1,gr_ngc
!!$ Vbeg(1:NUMB_VAR) = V(1:NUMB_VAR, gr_ibeg(XDIM)+i-1,j)
!!$ Vend(1:NUMB_VAR) = V(1:NUMB_VAR, gr_iend(XDIM)-i+1,j)
!!$
!!$ V(1:NUMB_VAR, gr_iend(XDIM)+i,j) = Vbeg(1:NUMB_VAR)
!!$ V(1:NUMB_VAR, gr_ibeg(XDIM)-i,j) = Vend(1:NUMB_VAR)
!!$ end do
!!$ end do
return
end subroutine bc_outflowX
subroutine bc_periodic(V)
!I need to think about how to deal with the corner cases only affects GP
implicit none
real, dimension(NUMB_VAR,gr_imax(XDIM),gr_imax(YDIM)), intent(INOUT) :: V
real, dimension(NUMB_VAR) :: Vbeg, Vend, Vne, Vnw, Vse, Vsw
integer :: i,j
!first do periodic in x
do j = gr_ibeg(YDIM),gr_iend(YDIM)!1,gr_ny
do i = 1,gr_ngc
Vbeg(1:NUMB_VAR) = V(1:NUMB_VAR, gr_ibeg(XDIM)+i-1,j)
Vend(1:NUMB_VAR) = V(1:NUMB_VAR, gr_iend(XDIM)-i+1,j)
V(1:NUMB_VAR, gr_iend(XDIM)+i,j) = Vbeg(1:NUMB_VAR)
V(1:NUMB_VAR, gr_ibeg(XDIM)-i,j) = Vend(1:NUMB_VAR)
end do
end do
!now periodic in y
do i = gr_ibeg(XDIM),gr_iend(XDIM)!1,gr_nx
do j = 1,gr_ngc
Vbeg(1:NUMB_VAR) = V(1:NUMB_VAR, i, gr_ibeg(YDIM)+j-1)
Vend(1:NUMB_VAR) = V(1:NUMB_VAR, i, gr_iend(YDIM) - j+1)
V(1:NUMB_VAR, i,gr_iend(YDIM)+j) = Vbeg(1:NUMB_VAR)
V(1:NUMB_VAR, i,gr_ibeg(YDIM)-j) = Vend(1:NUMB_VAR)
end do
end do
!deal with the corners
do i = 1, gr_ngc
do j = 1, gr_ngc
Vnw(:) = V(1:NUMB_VAR, gr_ibeg(XDIM) + i-1, gr_iend(YDIM) - j+1)
Vne(:) = V(1:NUMB_VAR, gr_iend(XDIM) - i+1, gr_iend(YDIM) -j+1)
Vsw(:) = V(1:NUMB_VAR, gr_ibeg(XDIM) + i-1, gr_ibeg(YDIM) + j-1)
Vse(:) = V(1:NUMB_VAR, gr_iend(XDIM) - i+1, gr_ibeg(YDIM) + j-1)
V(1:NUMB_VAR, gr_ibeg(XDIM) -i, gr_iend(YDIM) + j) = Vse !nw->se
V(1:NUMB_VAR, gr_iend(XDIM) +i, gr_iend(YDIM) + j) = Vsw !ne->sw
V(1:NUMB_VAR, gr_ibeg(XDIM) -i, gr_ibeg(YDIM) - j) = Vne !sw->ne
V(1:NUMB_VAR, gr_iend(XDIM) +i, gr_ibeg(YDIM) - j) = Vnw !se->nw
end do
end do
end subroutine bc_periodic
!!$ subroutine bc_user(V)
!!$ implicit none
!!$ real, dimension(NUMB_VAR,gr_imax), intent(INOUT) :: V
!!$ !BC for shu-osher problem.
!!$ !don't do anything!
!!$ end subroutine bc_user
!!$
!!$ subroutine bc_outflow(V)
!!$ implicit none
!!$ real, dimension(NUMB_VAR,gr_imax), intent(INOUT) :: V
!!$ integer :: i
!!$
!!$ do i = 1, gr_ngc
!!$ ! on the left GC
!!$ V(1:NUMB_VAR,i) = V(1:NUMB_VAR,gr_ibeg)
!!$
!!$ ! on the right GC
!!$ V(1:NUMB_VAR,gr_iend+i) = V(1:NUMB_VAR,gr_iend)
!!$ end do
!!$
!!$ return
!!$ end subroutine bc_outflow
!!$
!!$ subroutine bc_reflect(V)
!!$ implicit none
!!$ real, dimension(NUMB_VAR,gr_imax), intent(INOUT) :: V
!!$ integer :: i,k0,k1
!!$
!!$ do i = 1, gr_ngc
!!$ k0 = 2*gr_ngc+1
!!$ k1 = gr_iend-gr_ngc
!!$
!!$ ! on the left GC
!!$ V( :,i) = V( :,k0-i)
!!$ V(VELX_VAR,i) =-V(VELX_VAR,k0-i)
!!$
!!$ ! on the right GC
!!$ V( :,k1+k0-i) = V( :,k1+i)
!!$ V(VELX_VAR,k1+k0-i) =-V( VELX_VAR,k1+i)
!!$ end do
!!$
!!$ return
!!$ end subroutine bc_reflect
end module bc