MPI-AMRVAC  3.0
The MPI - Adaptive Mesh Refinement - Versatile Advection Code
mod_ghostcells_update.t
Go to the documentation of this file.
1 !> update ghost cells of all blocks including physical boundaries
3 
4  implicit none
5  ! Special buffer for pole copy
6  type wbuffer
7  double precision, dimension(:^D&,:), allocatable :: w
8  end type wbuffer
9 
10  ! A switch of update physical boundary or not
11  logical, public :: bcphys=.true.
12  ! A switch of update non-physical boundary or not
13  !logical, public :: bcexch=.true.
14  integer :: ixm^l, ixcog^l, ixcom^l
15 
16  ! The number of interleaving sending buffers for ghost cells
17  integer, parameter :: npwbuf=2
18 
19  ! The first index goes from -1:2, where -1 is used when a block touches the
20  ! lower boundary, 1 when a block touches an upper boundary, and 0 a situation
21  ! away from boundary conditions, 2 when a block touched both lower and upper
22  ! boundary
23 
24  ! index ranges to send (S) to sibling blocks, receive (R) from sibling blocks
25  integer, dimension(-1:2,-1:1) :: ixs_srl_^l, ixr_srl_^l
26 
27  ! index ranges of staggered variables to send (S) to sibling blocks, receive (R) from sibling blocks
28  integer, dimension(^ND,-1:1) :: ixs_srl_stg_^l, ixr_srl_stg_^l
29 
30  ! index ranges to send (S) restricted (r) ghost cells to coarser blocks
31  integer, dimension(-1:1,-1:1) :: ixs_r_^l
32 
33  ! index ranges of staggered variables to send (S) restricted (r) ghost cells to coarser blocks
34  integer, dimension(^ND,-1:1) :: ixs_r_stg_^l
35 
36  ! index ranges to receive restriced ghost cells from finer blocks
37  integer, dimension(-1:1, 0:3) :: ixr_r_^l
38 
39  ! index ranges of staggered variables to receive restriced ghost cells from finer blocks
40  integer, dimension(^ND,0:3) :: ixr_r_stg_^l
41 
42  ! send prolongated (p) ghost cells to finer blocks, receive prolongated from coarser blocks
43  integer, dimension(-1:1, 0:3) :: ixs_p_^l, ixr_p_^l
44 
45  ! send prolongated (p) staggered ghost cells to finer blocks, receive prolongated from coarser blocks
46  integer, dimension(^ND,0:3) :: ixs_p_stg_^l, ixr_p_stg_^l
47 
48  ! number of MPI receive-send pairs, srl: same refinement level; r: restrict; p: prolong
50 
51  ! record index position of buffer arrays
53 
54  ! count of times of send and receive
56 
57  ! count of times of send and receive for cell center ghost cells
58  integer :: isend_c, irecv_c
59 
60  ! tag of MPI send and recv
61  integer, private :: itag
62 
63  ! total sizes = cell-center normal flux + stagger-grid flux of send and receive
64  integer, dimension(-1:1^D&) :: sizes_srl_send_total, sizes_srl_recv_total
65 
66  ! sizes of buffer arrays for center-grid variable for siblings and restrict
67  integer, dimension(:), allocatable :: recvrequest_c_sr, sendrequest_c_sr
68  integer, dimension(:,:), allocatable :: recvstatus_c_sr, sendstatus_c_sr
69 
70  ! sizes of buffer arrays for center-grid variable for prolongation
71  integer, dimension(:), allocatable :: recvrequest_c_p, sendrequest_c_p
72  integer, dimension(:,:), allocatable :: recvstatus_c_p, sendstatus_c_p
73 
74  ! sizes of buffer arrays for stagger-grid variable
75  integer, dimension(^ND,-1:1^D&) :: sizes_srl_send_stg, sizes_srl_recv_stg
76 
77  integer, dimension(:), allocatable :: recvrequest_srl, sendrequest_srl
78  integer, dimension(:,:), allocatable :: recvstatus_srl, sendstatus_srl
79 
80  ! buffer arrays for send and receive of siblings, allocate in build_connectivity
81  double precision, dimension(:), allocatable :: recvbuffer_srl, sendbuffer_srl
82 
83  integer, dimension(:), allocatable :: recvrequest_r, sendrequest_r
84  integer, dimension(:,:), allocatable :: recvstatus_r, sendstatus_r
85 
86  ! buffer arrays for send and receive in restriction
87  double precision, dimension(:), allocatable :: recvbuffer_r, sendbuffer_r
88 
89  integer, dimension(:), allocatable :: recvrequest_p, sendrequest_p
90  integer, dimension(:,:), allocatable :: recvstatus_p, sendstatus_p
91 
92  ! buffer arrays for send and receive in prolongation
93  double precision, dimension(:), allocatable :: recvbuffer_p, sendbuffer_p
94 
95  ! sizes to allocate buffer arrays for send and receive for restriction
96  integer, dimension(-1:1^D&) :: sizes_r_send_total
97  integer, dimension(0:3^D&) :: sizes_r_recv_total
98  integer, dimension(^ND,-1:1^D&) :: sizes_r_send_stg
99  integer, dimension(^ND,0:3^D&) :: sizes_r_recv_stg
100 
101  ! sizes to allocate buffer arrays for send and receive for restriction
102  integer, dimension(0:3^D&) :: sizes_p_send_total, sizes_p_recv_total
103  integer, dimension(^ND,0:3^D&) :: sizes_p_send_stg, sizes_p_recv_stg
104 
105  ! There are two variants, _f indicates that all flux variables are filled,
106  ! whereas _p means that part of the variables is filled
107  ! Furthermore _r_ stands for restrict, _p_ for prolongation.
108  integer, dimension(-1:2^D&,-1:1^D&), target :: type_send_srl_f, type_recv_srl_f
109  integer, dimension(-1:1^D&,-1:1^D&), target :: type_send_r_f
110  integer, dimension(-1:1^D&, 0:3^D&), target :: type_recv_r_f, type_send_p_f, type_recv_p_f
111  integer, dimension(-1:2^D&,-1:1^D&), target :: type_send_srl_p1, type_recv_srl_p1
112  integer, dimension(-1:1^D&,-1:1^D&), target :: type_send_r_p1
113  integer, dimension(-1:1^D&, 0:3^D&), target :: type_recv_r_p1, type_send_p_p1, type_recv_p_p1
114  integer, dimension(-1:2^D&,-1:1^D&), target :: type_send_srl_p2, type_recv_srl_p2
115  integer, dimension(-1:1^D&,-1:1^D&), target :: type_send_r_p2
116  integer, dimension(-1:1^D&, 0:3^D&), target :: type_recv_r_p2, type_send_p_p2, type_recv_p_p2
117  integer, dimension(:^D&,:^D&), pointer :: type_send_srl, type_recv_srl, type_send_r
118  integer, dimension(:^D&,:^D&), pointer :: type_recv_r, type_send_p, type_recv_p
119 
120 contains
121 
122  subroutine init_bc()
125 
126  integer :: nghostcellsCo, interpolation_order
127  integer :: nx^D, nxCo^D, ixG^L, i^D, ic^D, inc^D, idir
128 
129  ixg^l=ixg^ll;
130  ixm^l=ixg^l^lsubnghostcells;
131  ixcogmin^d=1;
132  !ixCoGmax^D=ixGmax^D/2+nghostcells;
133  ixcogmax^d=(ixghi^d-2*nghostcells)/2+2*nghostcells;
134 
135  ixcom^l=ixcog^l^lsubnghostcells;
136 
137  nx^d=ixmmax^d-ixmmin^d+1;
138  nxco^d=nx^d/2;
139 
140  if(ghost_copy) then
141  interpolation_order=1
142  else
143  interpolation_order=2
144  end if
145  nghostcellsco=int((nghostcells+1)/2)
146 
147  if (nghostcellsco+interpolation_order-1>nghostcells) then
148  call mpistop("interpolation order for prolongation in getbc too high")
149  end if
150 
151  ! (iib,i) index has following meanings: iib = 0 means it is not at any physical boundary
152  ! iib=-1 means it is at the minimum side of a physical boundary
153  ! iib= 1 means it is at the maximum side of a physical boundary
154  ! i=-1 means subregion prepared for the neighbor at its minimum side
155  ! i= 1 means subregion prepared for the neighbor at its maximum side
156  {
157  ixs_srl_min^d(:,-1)=ixmmin^d
158  ixs_srl_min^d(:, 0)=ixmmin^d
159  ixs_srl_min^d(:, 1)=ixmmax^d+1-nghostcells
160  ixs_srl_max^d(:,-1)=ixmmin^d-1+nghostcells
161  ixs_srl_max^d(:, 0)=ixmmax^d
162  ixs_srl_max^d(:, 1)=ixmmax^d
163 
164  ixr_srl_min^d(:,-1)=1
165  ixr_srl_min^d(:, 0)=ixmmin^d
166  ixr_srl_min^d(:, 1)=ixmmax^d+1
167  ixr_srl_max^d(:,-1)=nghostcells
168  ixr_srl_max^d(:, 0)=ixmmax^d
169  ixr_srl_max^d(:, 1)=ixgmax^d
170 
171  ixs_r_min^d(:,-1)=ixcommin^d
172  ixs_r_min^d(:, 0)=ixcommin^d
173  ixs_r_min^d(:, 1)=ixcommax^d+1-nghostcells
174  ixs_r_max^d(:,-1)=ixcommin^d-1+nghostcells
175  ixs_r_max^d(:, 0)=ixcommax^d
176  ixs_r_max^d(:, 1)=ixcommax^d
177 
178  ixr_r_min^d(:, 0)=1
179  ixr_r_min^d(:, 1)=ixmmin^d
180  ixr_r_min^d(:, 2)=ixmmin^d+nxco^d
181  ixr_r_min^d(:, 3)=ixmmax^d+1
182  ixr_r_max^d(:, 0)=nghostcells
183  ixr_r_max^d(:, 1)=ixmmin^d-1+nxco^d
184  ixr_r_max^d(:, 2)=ixmmax^d
185  ixr_r_max^d(:, 3)=ixgmax^d
186 
187  ixs_p_min^d(:, 0)=ixmmin^d-(interpolation_order-1)
188  ixs_p_min^d(:, 1)=ixmmin^d-(interpolation_order-1)
189  ixs_p_min^d(:, 2)=ixmmin^d+nxco^d-nghostcellsco-(interpolation_order-1)
190  ixs_p_min^d(:, 3)=ixmmax^d+1-nghostcellsco-(interpolation_order-1)
191  ixs_p_max^d(:, 0)=ixmmin^d-1+nghostcellsco+(interpolation_order-1)
192  ixs_p_max^d(:, 1)=ixmmin^d-1+nxco^d+nghostcellsco+(interpolation_order-1)
193  ixs_p_max^d(:, 2)=ixmmax^d+(interpolation_order-1)
194  ixs_p_max^d(:, 3)=ixmmax^d+(interpolation_order-1)
195 
196  if(.not.phys_req_diagonal) then
197  ! exclude ghost-cell region when diagonal cells are unknown
198  ixs_p_min^d(:, 0)=ixmmin^d
199  ixs_p_max^d(:, 3)=ixmmax^d
200  ixs_p_max^d(:, 1)=ixmmin^d-1+nxco^d+(interpolation_order-1)
201  ixs_p_min^d(:, 2)=ixmmin^d+nxco^d-(interpolation_order-1)
202  end if
203 
204  ixr_p_min^d(:, 0)=ixcommin^d-nghostcellsco-(interpolation_order-1)
205  ixr_p_min^d(:, 1)=ixcommin^d-(interpolation_order-1)
206  ixr_p_min^d(:, 2)=ixcommin^d-nghostcellsco-(interpolation_order-1)
207  ixr_p_min^d(:, 3)=ixcommax^d+1-(interpolation_order-1)
208  ixr_p_max^d(:, 0)=nghostcells+(interpolation_order-1)
209  ixr_p_max^d(:, 1)=ixcommax^d+nghostcellsco+(interpolation_order-1)
210  ixr_p_max^d(:, 2)=ixcommax^d+(interpolation_order-1)
211  ixr_p_max^d(:, 3)=ixcommax^d+nghostcellsco+(interpolation_order-1)
212 
213  if(.not.phys_req_diagonal) then
214  ixr_p_max^d(:, 0)=nghostcells
215  ixr_p_min^d(:, 3)=ixcommax^d+1
216  ixr_p_max^d(:, 1)=ixcommax^d+(interpolation_order-1)
217  ixr_p_min^d(:, 2)=ixcommin^d-(interpolation_order-1)
218  end if
219 
220  \}
221 
222  if (stagger_grid) then
223  allocate(pole_buf%ws(ixgs^t,nws))
224  ! Staggered (face-allocated) variables
225  do idir=1,ndim
226  { ixs_srl_stg_min^d(idir,-1)=ixmmin^d-kr(idir,^d)
227  ixs_srl_stg_max^d(idir,-1)=ixmmin^d-1+nghostcells
228  ixs_srl_stg_min^d(idir,0) =ixmmin^d-kr(idir,^d)
229  ixs_srl_stg_max^d(idir,0) =ixmmax^d
230  ixs_srl_stg_min^d(idir,1) =ixmmax^d-nghostcells+1-kr(idir,^d)
231  ixs_srl_stg_max^d(idir,1) =ixmmax^d
232 
233  ixr_srl_stg_min^d(idir,-1)=1-kr(idir,^d)
234  ixr_srl_stg_max^d(idir,-1)=nghostcells
235  ixr_srl_stg_min^d(idir,0) =ixmmin^d-kr(idir,^d)
236  ixr_srl_stg_max^d(idir,0) =ixmmax^d
237  ixr_srl_stg_min^d(idir,1) =ixmmax^d+1-kr(idir,^d)
238  ixr_srl_stg_max^d(idir,1) =ixgmax^d
239 
240  ixs_r_stg_min^d(idir,-1)=ixcommin^d-kr(idir,^d)
241  ixs_r_stg_max^d(idir,-1)=ixcommin^d-1+nghostcells
242  ixs_r_stg_min^d(idir,0) =ixcommin^d-kr(idir,^d)
243  ixs_r_stg_max^d(idir,0) =ixcommax^d
244  ixs_r_stg_min^d(idir,1) =ixcommax^d+1-nghostcells-kr(idir,^d)
245  ixs_r_stg_max^d(idir,1) =ixcommax^d
246 
247  ixr_r_stg_min^d(idir,0)=1-kr(idir,^d)
248  ixr_r_stg_max^d(idir,0)=nghostcells
249  ixr_r_stg_min^d(idir,1)=ixmmin^d-kr(idir,^d)
250  ixr_r_stg_max^d(idir,1)=ixmmin^d-1+nxco^d
251  ixr_r_stg_min^d(idir,2)=ixmmin^d+nxco^d-kr(idir,^d)
252  ixr_r_stg_max^d(idir,2)=ixmmax^d
253  ixr_r_stg_min^d(idir,3)=ixmmax^d+1-kr(idir,^d)
254  ixr_r_stg_max^d(idir,3)=ixgmax^d
255  \}
256  {if (idir==^d) then
257  ! Parallel components
258  {
259  ixs_p_stg_min^d(idir,0)=ixmmin^d-1 ! -1 to make redundant
260  ixs_p_stg_max^d(idir,0)=ixmmin^d-1+nghostcellsco
261  ixs_p_stg_min^d(idir,1)=ixmmin^d-1 ! -1 to make redundant
262  ixs_p_stg_max^d(idir,1)=ixmmin^d-1+nxco^d+nghostcellsco
263  ixs_p_stg_min^d(idir,2)=ixmmax^d-nxco^d-nghostcellsco
264  ixs_p_stg_max^d(idir,2)=ixmmax^d
265  ixs_p_stg_min^d(idir,3)=ixmmax^d-nghostcellsco
266  ixs_p_stg_max^d(idir,3)=ixmmax^d
267 
268  ixr_p_stg_min^d(idir,0)=ixcommin^d-1-nghostcellsco
269  ixr_p_stg_max^d(idir,0)=ixcommin^d-1
270  ixr_p_stg_min^d(idir,1)=ixcommin^d-1 ! -1 to make redundant
271  ixr_p_stg_max^d(idir,1)=ixcommax^d+nghostcellsco
272  ixr_p_stg_min^d(idir,2)=ixcommin^d-1-nghostcellsco
273  ixr_p_stg_max^d(idir,2)=ixcommax^d
274  ixr_p_stg_min^d(idir,3)=ixcommax^d+1-1 ! -1 to make redundant
275  ixr_p_stg_max^d(idir,3)=ixcommax^d+nghostcellsco
276  \}
277  else
278  {
279  ! Perpendicular component
280  ixs_p_stg_min^d(idir,0)=ixmmin^d
281  ixs_p_stg_max^d(idir,0)=ixmmin^d-1+nghostcellsco+(interpolation_order-1)
282  ixs_p_stg_min^d(idir,1)=ixmmin^d
283  ixs_p_stg_max^d(idir,1)=ixmmin^d-1+nxco^d+nghostcellsco+(interpolation_order-1)
284  ixs_p_stg_min^d(idir,2)=ixmmax^d+1-nxco^d-nghostcellsco-(interpolation_order-1)
285  ixs_p_stg_max^d(idir,2)=ixmmax^d
286  ixs_p_stg_min^d(idir,3)=ixmmax^d+1-nghostcellsco-(interpolation_order-1)
287  ixs_p_stg_max^d(idir,3)=ixmmax^d
288 
289  ixr_p_stg_min^d(idir,0)=ixcommin^d-nghostcellsco-(interpolation_order-1)
290  ixr_p_stg_max^d(idir,0)=ixcommin^d-1
291  ixr_p_stg_min^d(idir,1)=ixcommin^d
292  ixr_p_stg_max^d(idir,1)=ixcommax^d+nghostcellsco+(interpolation_order-1)
293  ixr_p_stg_min^d(idir,2)=ixcommin^d-nghostcellsco-(interpolation_order-1)
294  ixr_p_stg_max^d(idir,2)=ixcommax^d
295  ixr_p_stg_min^d(idir,3)=ixcommax^d+1
296  ixr_p_stg_max^d(idir,3)=ixcommax^d+nghostcellsco+(interpolation_order-1)
297  \}
298  end if
299  }
300  end do
301  ! calculate sizes for buffer arrays for siblings
302  {do i^db=-1,1\}
303  ! Staggered (face-allocated) variables
304  do idir=1,ndim
305  sizes_srl_send_stg(idir,i^d)={(ixs_srl_stg_max^d(idir,i^d)-ixs_srl_stg_min^d(idir,i^d)+1)|*}
306  sizes_srl_recv_stg(idir,i^d)={(ixr_srl_stg_max^d(idir,i^d)-ixr_srl_stg_min^d(idir,i^d)+1)|*}
307  sizes_r_send_stg(idir,i^d)={(ixs_r_stg_max^d(idir,i^d)-ixs_r_stg_min^d(idir,i^d)+1)|*}
308  end do
311  sizes_r_send_total(i^d)=sum(sizes_r_send_stg(:,i^d))
312  {end do\}
313 
314  {do i^db=0,3\}
315  ! Staggered (face-allocated) variables
316  do idir=1,ndim
317  sizes_r_recv_stg(idir,i^d)={(ixr_r_stg_max^d(idir,i^d)-ixr_r_stg_min^d(idir,i^d)+1)|*}
318  sizes_p_send_stg(idir,i^d)={(ixs_p_stg_max^d(idir,i^d)-ixs_p_stg_min^d(idir,i^d)+1)|*}
319  sizes_p_recv_stg(idir,i^d)={(ixr_p_stg_max^d(idir,i^d)-ixr_p_stg_min^d(idir,i^d)+1)|*}
320  end do
321  sizes_r_recv_total(i^d)=sum(sizes_r_recv_stg(:,i^d))
322  sizes_p_send_total(i^d)=sum(sizes_p_send_stg(:,i^d))
323  sizes_p_recv_total(i^d)=sum(sizes_p_recv_stg(:,i^d))
324  {end do\}
325  end if
326  if(.not.stagger_grid .or. physics_type=='mf') then
327  ! extend index range to physical boundary
328  {
329  ixs_srl_min^d(-1,0)=1
330  ixs_srl_min^d( 1,0)=ixmmin^d
331  ixs_srl_min^d( 2,0)=1
332  ixs_srl_max^d(-1,0)=ixmmax^d
333  ixs_srl_max^d( 1,0)=ixgmax^d
334  ixs_srl_max^d( 2,0)=ixgmax^d
335 
336  ixr_srl_min^d(-1,0)=1
337  ixr_srl_min^d( 1,0)=ixmmin^d
338  ixr_srl_min^d( 2,0)=1
339  ixr_srl_max^d(-1,0)=ixmmax^d
340  ixr_srl_max^d( 1,0)=ixgmax^d
341  ixr_srl_max^d( 2,0)=ixgmax^d
342 
343  ixs_r_min^d(-1,0)=1
344  ixs_r_min^d( 1,0)=ixcommin^d
345  ixs_r_max^d(-1,0)=ixcommax^d
346  ixs_r_max^d( 1,0)=ixcogmax^d
347 
348  ixr_r_min^d(-1,1)=1
349  ixr_r_max^d(-1,1)=ixmmin^d-1+nxco^d
350  ixr_r_min^d( 1,2)=ixmmin^d+nxco^d
351  ixr_r_max^d( 1,2)=ixgmax^d
352 
353  ixs_p_min^d(-1,1)=1
354  ixs_p_max^d( 1,2)=ixgmax^d
355 
356  ixr_p_min^d(-1,1)=1
357  ixr_p_max^d( 1,2)=ixcogmax^d
358  \}
359  end if
360 
361  end subroutine init_bc
362 
363  subroutine create_bc_mpi_datatype(nwstart,nwbc)
365 
366  integer, intent(in) :: nwstart, nwbc
367  integer :: i^D, ic^D, inc^D, iib^D
368 
369  {do i^db=-1,1\}
370  if (i^d==0|.and.) cycle
371  {do iib^db=-1,2\}
372  call get_bc_comm_type(type_send_srl(iib^d,i^d),ixs_srl_^l(iib^d,i^d),ixg^ll,nwstart,nwbc)
373  call get_bc_comm_type(type_recv_srl(iib^d,i^d),ixr_srl_^l(iib^d,i^d),ixg^ll,nwstart,nwbc)
374  if (iib^d==2|.or.) cycle
375  call get_bc_comm_type(type_send_r(iib^d,i^d), ixs_r_^l(iib^d,i^d),ixcog^l,nwstart,nwbc)
376  {do ic^db=1+int((1-i^db)/2),2-int((1+i^db)/2)
377  inc^db=2*i^db+ic^db\}
378  call get_bc_comm_type(type_recv_r(iib^d,inc^d),ixr_r_^l(iib^d,inc^d), ixg^ll,nwstart,nwbc)
379  call get_bc_comm_type(type_send_p(iib^d,inc^d),ixs_p_^l(iib^d,inc^d), ixg^ll,nwstart,nwbc)
380  call get_bc_comm_type(type_recv_p(iib^d,inc^d),ixr_p_^l(iib^d,inc^d),ixcog^l,nwstart,nwbc)
381  {end do\}
382  {end do\}
383  {end do\}
384 
385  end subroutine create_bc_mpi_datatype
386 
387  subroutine get_bc_comm_type(comm_type,ix^L,ixG^L,nwstart,nwbc)
389 
390  integer, intent(inout) :: comm_type
391  integer, intent(in) :: ix^L, ixG^L, nwstart, nwbc
392 
393  integer, dimension(ndim+1) :: fullsize, subsize, start
394 
395  ^d&fullsize(^d)=ixgmax^d;
396  fullsize(ndim+1)=nw
397  ^d&subsize(^d)=ixmax^d-ixmin^d+1;
398  subsize(ndim+1)=nwbc
399  ^d&start(^d)=ixmin^d-1;
400  start(ndim+1)=nwstart-1
401 
402  call mpi_type_create_subarray(ndim+1,fullsize,subsize,start,mpi_order_fortran, &
403  mpi_double_precision,comm_type,ierrmpi)
404  call mpi_type_commit(comm_type,ierrmpi)
405 
406  end subroutine get_bc_comm_type
407 
408  subroutine put_bc_comm_types()
410 
411  integer :: i^D, ic^D, inc^D, iib^D
412 
413  {do i^db=-1,1\}
414  if (i^d==0|.and.) cycle
415  {do iib^db=-1,2\}
416  call mpi_type_free(type_send_srl(iib^d,i^d),ierrmpi)
417  call mpi_type_free(type_recv_srl(iib^d,i^d),ierrmpi)
418  if (levmin==levmax) cycle
419  if (iib^d==2|.or.) cycle
420  call mpi_type_free(type_send_r(iib^d,i^d),ierrmpi)
421  {do ic^db=1+int((1-i^db)/2),2-int((1+i^db)/2)
422  inc^db=2*i^db+ic^db\}
423  call mpi_type_free(type_recv_r(iib^d,inc^d),ierrmpi)
424  call mpi_type_free(type_send_p(iib^d,inc^d),ierrmpi)
425  call mpi_type_free(type_recv_p(iib^d,inc^d),ierrmpi)
426  {end do\}
427  {end do\}
428  {end do\}
429 
430  end subroutine put_bc_comm_types
431 
432  !> do update ghost cells of all blocks including physical boundaries
433  subroutine getbc(time,qdt,psb,nwstart,nwbc,req_diag)
435  use mod_physics
436 
437  double precision, intent(in) :: time, qdt
438  type(state), target :: psb(max_blocks)
439  integer, intent(in) :: nwstart ! Fill from nwstart
440  integer, intent(in) :: nwbc ! Number of variables to fill
441  logical, intent(in), optional :: req_diag ! If false, skip diagonal ghost cells
442 
443  double precision :: time_bcin
444  integer :: ipole, nwhead, nwtail
445  integer :: iigrid, igrid, ineighbor, ipe_neighbor, isizes
446  integer :: ixR^L, ixS^L
447  integer :: i^D, n_i^D, ic^D, inc^D, n_inc^D, iib^D, idir
448  ! store physical boundary indicating index
449  integer :: idphyb(ndim,max_blocks)
450  integer :: isend_buf(npwbuf), ipwbuf, nghostcellsco
451  ! index pointer for buffer arrays as a start for a segment
452  integer :: ibuf_start, ibuf_next
453  ! shapes of reshape
454  integer, dimension(1) :: shapes
455  logical :: req_diagonal
456  type(wbuffer) :: pwbuf(npwbuf)
457 
458  time_bcin=mpi_wtime()
459 
460  nwhead=nwstart
461  nwtail=nwstart+nwbc-1
462 
463  req_diagonal = .true.
464  if (present(req_diag)) req_diagonal = req_diag
465 
466  ! fill internal physical boundary
467  if (internalboundary) then
468  call getintbc(time,ixg^ll)
469  end if
470 
471  ! fill physical-boundary ghost cells before internal ghost-cell values exchange
472  if(bcphys.and. .not.stagger_grid) then
473  !$OMP PARALLEL DO SCHEDULE(dynamic) PRIVATE(igrid)
474  do iigrid=1,igridstail; igrid=igrids(iigrid);
475  if(.not.phyboundblock(igrid)) cycle
476  call fill_boundary_before_gc(igrid)
477  end do
478  !$OMP END PARALLEL DO
479  end if
480 
481  ! prepare coarse values to send to coarser neighbors
482  !$OMP PARALLEL DO SCHEDULE(dynamic) PRIVATE(igrid)
483  do iigrid=1,igridstail; igrid=igrids(iigrid);
484  if(any(neighbor_type(:^d&,igrid)==neighbor_coarse)) then
485  call coarsen_grid(psb(igrid),ixg^ll,ixm^l,psc(igrid),ixcog^l,ixcom^l)
486  {do i^db=-1,1\}
487  if(skip_direction([ i^d ])) cycle
488  if(neighbor_type(i^d,igrid)==neighbor_coarse) call fill_coarse_boundary(igrid,i^d)
489  {end do\}
490  end if
491  end do
492  !$OMP END PARALLEL DO
493  !if(bcexch) then
494  ! default : no singular axis
495  ipole=0
496  irecv_c=0
497  isend_c=0
498  isend_buf=0
499  ipwbuf=1
500 
501  if(stagger_grid) then
502  ibuf_recv_srl=1
503  ibuf_recv_r=1
504  ibuf_recv_p=1
505  ibuf_send_srl=1
506  ibuf_send_r=1
507  ibuf_send_p=1
508  irecv_srl=0
509  irecv_r=0
510  irecv_p=0
511  isend_srl=0
512  isend_r=0
513  isend_p=0
514  end if
515 
516 
517  ! MPI receive ghost-cell values from sibling blocks and finer neighbors in different processors
518  do iigrid=1,igridstail; igrid=igrids(iigrid);
519  call identifyphysbound(ps(igrid),iib^d)
520  ^d&idphyb(^d,igrid)=iib^d;
521  {do i^db=-1,1\}
522  if (skip_direction([ i^d ])) cycle
523  select case (neighbor_type(i^d,igrid))
524  case (neighbor_sibling)
525  call bc_recv_srl
526  case (neighbor_fine)
527  call bc_recv_restrict
528  end select
529  {end do\}
530  end do
531 
532  ! MPI send ghost-cell values to sibling blocks and coarser neighbors in different processors
533  do iigrid=1,igridstail; igrid=igrids(iigrid);
534  ^d&iib^d=idphyb(^d,igrid);
535  {do i^db=-1,1\}
536  if(skip_direction([ i^d ])) cycle
537  select case (neighbor_type(i^d,igrid))
538  case (neighbor_sibling)
539  call bc_send_srl
540  case (neighbor_coarse)
541  call bc_send_restrict
542  end select
543  {end do\}
544  end do
545 
546  ! fill ghost-cell values of sibling blocks and coarser neighbors in the same processor
547  !$OMP PARALLEL DO SCHEDULE(dynamic) PRIVATE(igrid,iib^D)
548  do iigrid=1,igridstail; igrid=igrids(iigrid);
549  ^d&iib^d=idphyb(^d,igrid);
550  {do i^db=-1,1\}
551  if(skip_direction([ i^d ])) cycle
552  select case (neighbor_type(i^d,igrid))
553  case(neighbor_sibling)
554  call bc_fill_srl(igrid,i^d,iib^d)
555  case(neighbor_coarse)
556  call bc_fill_restrict(igrid,i^d,iib^d)
557  end select
558  {end do\}
559  end do
560  !$OMP END PARALLEL DO
561 
562  call mpi_waitall(irecv_c,recvrequest_c_sr,recvstatus_c_sr,ierrmpi)
563  call mpi_waitall(isend_c,sendrequest_c_sr,sendstatus_c_sr,ierrmpi)
564 
565  if(stagger_grid) then
566  call mpi_waitall(nrecv_bc_srl,recvrequest_srl,recvstatus_srl,ierrmpi)
567  call mpi_waitall(nsend_bc_srl,sendrequest_srl,sendstatus_srl,ierrmpi)
568  call mpi_waitall(nrecv_bc_r,recvrequest_r,recvstatus_r,ierrmpi)
569  call mpi_waitall(nsend_bc_r,sendrequest_r,sendstatus_r,ierrmpi)
570  ! unpack the received data from sibling blocks and finer neighbors to fill ghost-cell staggered values
571  ibuf_recv_srl=1
572  ibuf_recv_r=1
573  do iigrid=1,igridstail; igrid=igrids(iigrid);
574  ^d&iib^d=idphyb(^d,igrid);
575  {do i^db=-1,1\}
576  if (skip_direction([ i^d ])) cycle
577  select case (neighbor_type(i^d,igrid))
578  case (neighbor_sibling)
579  call bc_fill_srl_stg
580  case (neighbor_fine)
582  end select
583  {end do\}
584  end do
585  end if
586 
587  do ipwbuf=1,npwbuf
588  if (isend_buf(ipwbuf)/=0) deallocate(pwbuf(ipwbuf)%w)
589  end do
590 
591  irecv_c=0
592  isend_c=0
593  isend_buf=0
594  ipwbuf=1
595 
596  ! MPI receive ghost-cell values from coarser neighbors in different processors
597  do iigrid=1,igridstail; igrid=igrids(iigrid);
598  ^d&iib^d=idphyb(^d,igrid);
599  {do i^db=-1,1\}
600  if (skip_direction([ i^d ])) cycle
601  if (neighbor_type(i^d,igrid)==neighbor_coarse) call bc_recv_prolong
602  {end do\}
603  end do
604  ! MPI send ghost-cell values to finer neighbors in different processors
605  do iigrid=1,igridstail; igrid=igrids(iigrid);
606  ^d&iib^d=idphyb(^d,igrid);
607  {do i^db=-1,1\}
608  if (skip_direction([ i^d ])) cycle
609  if (neighbor_type(i^d,igrid)==neighbor_fine) call bc_send_prolong
610  {end do\}
611  end do
612 
613  ! fill coarse ghost-cell values of finer neighbors in the same processor
614  !$OMP PARALLEL DO SCHEDULE(dynamic) PRIVATE(igrid,iib^D)
615  do iigrid=1,igridstail; igrid=igrids(iigrid);
616  ^d&iib^d=idphyb(^d,igrid);
617  {do i^db=-1,1\}
618  if (skip_direction([ i^d ])) cycle
619  if (neighbor_type(i^d,igrid)==neighbor_fine) call bc_fill_prolong(igrid,i^d,iib^d)
620  {end do\}
621  end do
622  !$OMP END PARALLEL DO
623 
624  call mpi_waitall(irecv_c,recvrequest_c_p,recvstatus_c_p,ierrmpi)
625  call mpi_waitall(isend_c,sendrequest_c_p,sendstatus_c_p,ierrmpi)
626 
627  if(stagger_grid) then
628  call mpi_waitall(nrecv_bc_p,recvrequest_p,recvstatus_p,ierrmpi)
629  call mpi_waitall(nsend_bc_p,sendrequest_p,sendstatus_p,ierrmpi)
630 
631  ! fill coarser representative ghost cells after receipt
632  ibuf_recv_p=1
633  do iigrid=1,igridstail; igrid=igrids(iigrid);
634  ^d&iib^d=idphyb(^d,igrid);
635  {do i^db=-1,1\}
636  if (skip_direction([ i^d ])) cycle
637  if(neighbor_type(i^d,igrid)==neighbor_coarse) call bc_fill_prolong_stg
638  {end do\}
639  end do
640  end if
641  ! do prolongation on the ghost-cell values based on the received coarse values from coarser neighbors
642  !$OMP PARALLEL DO SCHEDULE(dynamic) PRIVATE(igrid)
643  do iigrid=1,igridstail; igrid=igrids(iigrid);
644  call gc_prolong(igrid)
645  end do
646  !$OMP END PARALLEL DO
647 
648  do ipwbuf=1,npwbuf
649  if (isend_buf(ipwbuf)/=0) deallocate(pwbuf(ipwbuf)%w)
650  end do
651  !endif !bcexch
652  ! fill physical boundary ghost cells after internal ghost-cell values exchange
653  if(bcphys.and.stagger_grid) then
654  !$OMP PARALLEL DO SCHEDULE(dynamic) PRIVATE(igrid)
655  do iigrid=1,igridstail; igrid=igrids(iigrid);
656  if(.not.phyboundblock(igrid)) cycle
657  call fill_boundary_after_gc(igrid)
658  end do
659  !$OMP END PARALLEL DO
660  end if
661 
662  ! modify normal component of magnetic field to fix divB=0
663  if(bcphys.and.associated(phys_boundary_adjust)) then
664  !$OMP PARALLEL DO SCHEDULE(dynamic) PRIVATE(igrid)
665  do iigrid=1,igridstail; igrid=igrids(iigrid);
666  if(.not.phyboundblock(igrid)) cycle
667  call phys_boundary_adjust(igrid,psb)
668  end do
669  !$OMP END PARALLEL DO
670  end if
671 
672  time_bc=time_bc+(mpi_wtime()-time_bcin)
673 
674  contains
675 
676  logical function skip_direction(dir)
677  integer, intent(in) :: dir(^nd)
678 
679  if (all(dir == 0)) then
680  skip_direction = .true.
681  else if (.not. req_diagonal .and. count(dir /= 0) > 1) then
682  skip_direction = .true.
683  else
684  skip_direction = .false.
685  end if
686  end function skip_direction
687 
688  !> Physical boundary conditions
689  subroutine fill_boundary_before_gc(igrid)
690 
691  integer, intent(in) :: igrid
692 
693  integer :: idims,iside,i^D,k^L,ixB^L
694 
695  block=>psb(igrid)
696  ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
697  do idims=1,ndim
698  ! to avoid using as yet unknown corner info in more than 1D, we
699  ! fill only interior mesh ranges of the ghost cell ranges at first,
700  ! and progressively enlarge the ranges to include corners later
701  {
702  kmin^d=merge(0, 1, idims==^d)
703  kmax^d=merge(0, 1, idims==^d)
704  ixbmin^d=ixglo^d+kmin^d*nghostcells
705  ixbmax^d=ixghi^d-kmax^d*nghostcells
706  \}
707  {^iftwod
708  if(idims > 1 .and. neighbor_type(-1,0,igrid)==neighbor_boundary) ixbmin1=ixglo1
709  if(idims > 1 .and. neighbor_type( 1,0,igrid)==neighbor_boundary) ixbmax1=ixghi1}
710  {^ifthreed
711  if(idims > 1 .and. neighbor_type(-1,0,0,igrid)==neighbor_boundary) ixbmin1=ixglo1
712  if(idims > 1 .and. neighbor_type( 1,0,0,igrid)==neighbor_boundary) ixbmax1=ixghi1
713  if(idims > 2 .and. neighbor_type(0,-1,0,igrid)==neighbor_boundary) ixbmin2=ixglo2
714  if(idims > 2 .and. neighbor_type(0, 1,0,igrid)==neighbor_boundary) ixbmax2=ixghi2}
715  do iside=1,2
716  i^d=kr(^d,idims)*(2*iside-3);
717  if (aperiodb(idims)) then
718  if (neighbor_type(i^d,igrid) /= neighbor_boundary .and. &
719  .not. psb(igrid)%is_physical_boundary(2*idims-2+iside)) cycle
720  else
721  if (neighbor_type(i^d,igrid) /= neighbor_boundary) cycle
722  end if
723  call bc_phys(iside,idims,time,qdt,psb(igrid),ixg^ll,ixb^l)
724  end do
725  end do
726 
727  end subroutine fill_boundary_before_gc
728 
729  !> Physical boundary conditions
730  subroutine fill_boundary_after_gc(igrid)
731 
732  integer, intent(in) :: igrid
733 
734  integer :: idims,iside,i^D,k^L,ixB^L
735 
736  block=>psb(igrid)
737  ^d&dxlevel(^d)=rnode(rpdx^d_,igrid);
738  do idims=1,ndim
739  ! to avoid using as yet unknown corner info in more than 1D, we
740  ! fill only interior mesh ranges of the ghost cell ranges at first,
741  ! and progressively enlarge the ranges to include corners later
742  kmin1=0; kmax1=0;
743  {^iftwod
744  kmin2=merge(1, 0, idims .lt. 2 .and. neighbor_type(0,-1,igrid)==1)
745  kmax2=merge(1, 0, idims .lt. 2 .and. neighbor_type(0, 1,igrid)==1)}
746  {^ifthreed
747  kmin2=merge(1, 0, idims .lt. 2 .and. neighbor_type(0,-1,0,igrid)==1)
748  kmax2=merge(1, 0, idims .lt. 2 .and. neighbor_type(0, 1,0,igrid)==1)
749  kmin3=merge(1, 0, idims .lt. 3 .and. neighbor_type(0,0,-1,igrid)==1)
750  kmax3=merge(1, 0, idims .lt. 3 .and. neighbor_type(0,0, 1,igrid)==1)}
751  ixbmin^d=ixglo^d+kmin^d*nghostcells;
752  ixbmax^d=ixghi^d-kmax^d*nghostcells;
753  do iside=1,2
754  i^d=kr(^d,idims)*(2*iside-3);
755  if (aperiodb(idims)) then
756  if (neighbor_type(i^d,igrid) /= neighbor_boundary .and. &
757  .not. psb(igrid)%is_physical_boundary(2*idims-2+iside)) cycle
758  else
759  if (neighbor_type(i^d,igrid) /= neighbor_boundary) cycle
760  end if
761  call bc_phys(iside,idims,time,qdt,psb(igrid),ixg^ll,ixb^l)
762  end do
763  end do
764 
765  end subroutine fill_boundary_after_gc
766 
767  !> Receive from sibling at same refinement level
768  subroutine bc_recv_srl
769 
770  ipe_neighbor=neighbor(2,i^d,igrid)
771  if (ipe_neighbor/=mype) then
772  irecv_c=irecv_c+1
773  itag=(3**^nd+4**^nd)*(igrid-1)+{(i^d+1)*3**(^d-1)+}
774  call mpi_irecv(psb(igrid)%w,1,type_recv_srl(iib^d,i^d), &
775  ipe_neighbor,itag,icomm,recvrequest_c_sr(irecv_c),ierrmpi)
776  if(stagger_grid) then
778  call mpi_irecv(recvbuffer_srl(ibuf_recv_srl),sizes_srl_recv_total(i^d),mpi_double_precision, &
779  ipe_neighbor,itag,icomm,recvrequest_srl(irecv_srl),ierrmpi)
781  end if
782  end if
783 
784  end subroutine bc_recv_srl
785 
786  !> Receive from fine neighbor
787  subroutine bc_recv_restrict
788 
789  {do ic^db=1+int((1-i^db)/2),2-int((1+i^db)/2)
790  inc^db=2*i^db+ic^db\}
791  ipe_neighbor=neighbor_child(2,inc^d,igrid)
792  if (ipe_neighbor/=mype) then
793  irecv_c=irecv_c+1
794  itag=(3**^nd+4**^nd)*(igrid-1)+3**^nd+{inc^d*4**(^d-1)+}
795  call mpi_irecv(psb(igrid)%w,1,type_recv_r(iib^d,inc^d), &
796  ipe_neighbor,itag,icomm,recvrequest_c_sr(irecv_c),ierrmpi)
797  if(stagger_grid) then
798  irecv_r=irecv_r+1
799  call mpi_irecv(recvbuffer_r(ibuf_recv_r),sizes_r_recv_total(inc^d), &
800  mpi_double_precision,ipe_neighbor,itag, &
801  icomm,recvrequest_r(irecv_r),ierrmpi)
803  end if
804  end if
805  {end do\}
806 
807  end subroutine bc_recv_restrict
808 
809  !> Send to sibling at same refinement level
810  subroutine bc_send_srl
811 
812  ipe_neighbor=neighbor(2,i^d,igrid)
813 
814  if(ipe_neighbor/=mype) then
815  ineighbor=neighbor(1,i^d,igrid)
816  ipole=neighbor_pole(i^d,igrid)
817  if(ipole==0) then
818  n_i^d=-i^d;
819  isend_c=isend_c+1
820  itag=(3**^nd+4**^nd)*(ineighbor-1)+{(n_i^d+1)*3**(^d-1)+}
821  call mpi_isend(psb(igrid)%w,1,type_send_srl(iib^d,i^d), &
822  ipe_neighbor,itag,icomm,sendrequest_c_sr(isend_c),ierrmpi)
823  if(stagger_grid) then
824  ibuf_start=ibuf_send_srl
825  do idir=1,ndim
826  ixs^l=ixs_srl_stg_^l(idir,i^d);
827  ibuf_next=ibuf_start+sizes_srl_send_stg(idir,i^d)
828  shapes=(/sizes_srl_send_stg(idir,i^d)/)
829  sendbuffer_srl(ibuf_start:ibuf_next-1)=&
830  reshape(psb(igrid)%ws(ixs^s,idir),shapes)
831  ibuf_start=ibuf_next
832  end do
834  call mpi_isend(sendbuffer_srl(ibuf_send_srl),sizes_srl_send_total(i^d),&
835  mpi_double_precision, ipe_neighbor,itag,icomm,sendrequest_srl(isend_srl),ierrmpi)
836  ibuf_send_srl=ibuf_next
837  end if
838  else
839  ixs^l=ixs_srl_^l(iib^d,i^d);
840  select case (ipole)
841  {case (^d)
842  n_i^d=i^d^d%n_i^dd=-i^dd;\}
843  end select
844  if (isend_buf(ipwbuf)/=0) then
845  call mpi_wait(sendrequest_c_sr(isend_buf(ipwbuf)), &
846  sendstatus_c_sr(:,isend_buf(ipwbuf)),ierrmpi)
847  deallocate(pwbuf(ipwbuf)%w)
848  end if
849  allocate(pwbuf(ipwbuf)%w(ixs^s,nwhead:nwtail))
850  call pole_buffer(pwbuf(ipwbuf)%w,ixs^l,ixs^l,psb(igrid)%w,ixg^ll,ixs^l)
851  isend_c=isend_c+1
852  isend_buf(ipwbuf)=isend_c
853  itag=(3**^nd+4**^nd)*(ineighbor-1)+{(n_i^d+1)*3**(^d-1)+}
854  isizes={(ixsmax^d-ixsmin^d+1)*}*nwbc
855  call mpi_isend(pwbuf(ipwbuf)%w,isizes,mpi_double_precision, &
856  ipe_neighbor,itag,icomm,sendrequest_c_sr(isend_c),ierrmpi)
857  ipwbuf=1+modulo(ipwbuf,npwbuf)
858  if(stagger_grid) then
859  ibuf_start=ibuf_send_srl
860  do idir=1,ndim
861  ixs^l=ixs_srl_stg_^l(idir,i^d);
862  ibuf_next=ibuf_start+sizes_srl_send_stg(idir,i^d)
863  shapes=(/sizes_srl_send_stg(idir,i^d)/)
864  sendbuffer_srl(ibuf_start:ibuf_next-1)=&
865  reshape(psb(igrid)%ws(ixs^s,idir),shapes)
866  ibuf_start=ibuf_next
867  end do
869  call mpi_isend(sendbuffer_srl(ibuf_send_srl),sizes_srl_send_total(i^d),&
870  mpi_double_precision, ipe_neighbor,itag,icomm,sendrequest_srl(isend_srl),ierrmpi)
871  ibuf_send_srl=ibuf_next
872  end if
873  end if
874  end if
875 
876  end subroutine bc_send_srl
877 
878  subroutine bc_fill_srl(igrid,i^D,iib^D)
879  integer, intent(in) :: igrid,i^D,iib^D
880  integer :: ineighbor,ipe_neighbor,ipole,ixS^L,ixR^L,n_i^D,idir
881 
882  ipe_neighbor=neighbor(2,i^d,igrid)
883  if(ipe_neighbor==mype) then
884  ineighbor=neighbor(1,i^d,igrid)
885  ipole=neighbor_pole(i^d,igrid)
886  if(ipole==0) then
887  n_i^d=-i^d;
888  ixs^l=ixs_srl_^l(iib^d,i^d);
889  ixr^l=ixr_srl_^l(iib^d,n_i^d);
890  psb(ineighbor)%w(ixr^s,nwhead:nwtail)=&
891  psb(igrid)%w(ixs^s,nwhead:nwtail)
892  if(stagger_grid) then
893  do idir=1,ndim
894  ixs^l=ixs_srl_stg_^l(idir,i^d);
895  ixr^l=ixr_srl_stg_^l(idir,n_i^d);
896  psb(ineighbor)%ws(ixr^s,idir)=psb(igrid)%ws(ixs^s,idir)
897  end do
898  end if
899  else
900  ixs^l=ixs_srl_^l(iib^d,i^d);
901  select case (ipole)
902  {case (^d)
903  n_i^d=i^d^d%n_i^dd=-i^dd;\}
904  end select
905  ixr^l=ixr_srl_^l(iib^d,n_i^d);
906  call pole_copy(psb(ineighbor)%w,ixg^ll,ixr^l,psb(igrid)%w,ixg^ll,ixs^l,ipole)
907  if(stagger_grid) then
908  do idir=1,ndim
909  ixs^l=ixs_srl_stg_^l(idir,i^d);
910  ixr^l=ixr_srl_stg_^l(idir,n_i^d);
911  call pole_copy_stg(psb(ineighbor)%ws,ixr^l,psb(igrid)%ws,ixs^l,idir,ipole)
912  end do
913  end if
914  end if
915  end if
916 
917  end subroutine bc_fill_srl
918 
919  subroutine fill_coarse_boundary(igrid,i^D)
920  integer, intent(in) :: igrid,i^D
921 
922  integer :: idims,iside,k^L,ixB^L,ii^D
923 
924  if(phyboundblock(igrid).and..not.stagger_grid.and.bcphys) then
925  ! to use block in physical boundary setup for coarse representative
926  block=>psc(igrid)
927  ! filling physical boundary ghost cells of a coarser representative block for
928  ! sending swap region with width of nghostcells to its coarser neighbor
929  do idims=1,ndim
930  ! to avoid using as yet unknown corner info in more than 1D, we
931  ! fill only interior mesh ranges of the ghost cell ranges at first,
932  ! and progressively enlarge the ranges to include corners later
933  {kmin^d=merge(0, 1, idims==^d)
934  kmax^d=merge(0, 1, idims==^d)
935  ixbmin^d=ixcogmin^d+kmin^d*nghostcells
936  ixbmax^d=ixcogmax^d-kmax^d*nghostcells\}
937  {^iftwod
938  if(idims > 1 .and. neighbor_type(-1,0,igrid)==neighbor_boundary) ixbmin1=ixcogmin1
939  if(idims > 1 .and. neighbor_type( 1,0,igrid)==neighbor_boundary) ixbmax1=ixcogmax1}
940  {^ifthreed
941  if(idims > 1 .and. neighbor_type(-1,0,0,igrid)==neighbor_boundary) ixbmin1=ixcogmin1
942  if(idims > 1 .and. neighbor_type( 1,0,0,igrid)==neighbor_boundary) ixbmax1=ixcogmax1
943  if(idims > 2 .and. neighbor_type(0,-1,0,igrid)==neighbor_boundary) ixbmin2=ixcogmin2
944  if(idims > 2 .and. neighbor_type(0, 1,0,igrid)==neighbor_boundary) ixbmax2=ixcogmax2}
945  {if(i^d==-1) then
946  ixbmin^d=ixcogmin^d+nghostcells
947  ixbmax^d=ixcogmin^d+2*nghostcells-1
948  else if(i^d==1) then
949  ixbmin^d=ixcogmax^d-2*nghostcells+1
950  ixbmax^d=ixcogmax^d-nghostcells
951  end if\}
952  do iside=1,2
953  ii^d=kr(^d,idims)*(2*iside-3);
954  if ({abs(i^d)==1.and.abs(ii^d)==1|.or.}) cycle
955  if (neighbor_type(ii^d,igrid)/=neighbor_boundary) cycle
956  call bc_phys(iside,idims,time,0.d0,psc(igrid),ixcog^l,ixb^l)
957  end do
958  end do
959  end if
960 
961  end subroutine fill_coarse_boundary
962 
963  !> Send to coarser neighbor
964  subroutine bc_send_restrict
965 
966  ipe_neighbor=neighbor(2,i^d,igrid)
967  if(ipe_neighbor/=mype) then
968  ic^d=1+modulo(node(pig^d_,igrid)-1,2);
969  if({.not.(i^d==0.or.i^d==2*ic^d-3)|.or.}) return
970  ineighbor=neighbor(1,i^d,igrid)
971  ipole=neighbor_pole(i^d,igrid)
972  if(ipole==0) then
973  n_inc^d=-2*i^d+ic^d;
974  isend_c=isend_c+1
975  itag=(3**^nd+4**^nd)*(ineighbor-1)+3**^nd+{n_inc^d*4**(^d-1)+}
976  call mpi_isend(psc(igrid)%w,1,type_send_r(iib^d,i^d), &
977  ipe_neighbor,itag,icomm,sendrequest_c_sr(isend_c),ierrmpi)
978  if(stagger_grid) then
979  ibuf_start=ibuf_send_r
980  do idir=1,ndim
981  ixs^l=ixs_r_stg_^l(idir,i^d);
982  ibuf_next=ibuf_start+sizes_r_send_stg(idir,i^d)
983  shapes=(/sizes_r_send_stg(idir,i^d)/)
984  sendbuffer_r(ibuf_start:ibuf_next-1)=&
985  reshape(psc(igrid)%ws(ixs^s,idir),shapes)
986  ibuf_start=ibuf_next
987  end do
988  isend_r=isend_r+1
989  call mpi_isend(sendbuffer_r(ibuf_send_r),sizes_r_send_total(i^d),&
990  mpi_double_precision,ipe_neighbor,itag, &
991  icomm,sendrequest_r(isend_r),ierrmpi)
992  ibuf_send_r=ibuf_next
993  end if
994  else
995  ixs^l=ixs_r_^l(iib^d,i^d);
996  select case (ipole)
997  {case (^d)
998  n_inc^d=2*i^d+(3-ic^d)^d%n_inc^dd=-2*i^dd+ic^dd;\}
999  end select
1000  if(isend_buf(ipwbuf)/=0) then
1001  call mpi_wait(sendrequest_c_sr(isend_buf(ipwbuf)), &
1002  sendstatus_c_sr(:,isend_buf(ipwbuf)),ierrmpi)
1003  deallocate(pwbuf(ipwbuf)%w)
1004  end if
1005  allocate(pwbuf(ipwbuf)%w(ixs^s,nwhead:nwtail))
1006  call pole_buffer(pwbuf(ipwbuf)%w,ixs^l,ixs^l,psc(igrid)%w,ixcog^l,ixs^l)
1007  isend_c=isend_c+1
1008  isend_buf(ipwbuf)=isend_c
1009  itag=(3**^nd+4**^nd)*(ineighbor-1)+3**^nd+{n_inc^d*4**(^d-1)+}
1010  isizes={(ixsmax^d-ixsmin^d+1)*}*nwbc
1011  call mpi_isend(pwbuf(ipwbuf)%w,isizes,mpi_double_precision, &
1012  ipe_neighbor,itag,icomm,sendrequest_c_sr(isend_c),ierrmpi)
1013  ipwbuf=1+modulo(ipwbuf,npwbuf)
1014  if(stagger_grid) then
1015  ibuf_start=ibuf_send_r
1016  do idir=1,ndim
1017  ixs^l=ixs_r_stg_^l(idir,i^d);
1018  ibuf_next=ibuf_start+sizes_r_send_stg(idir,i^d)
1019  shapes=(/sizes_r_send_stg(idir,i^d)/)
1020  sendbuffer_r(ibuf_start:ibuf_next-1)=&
1021  reshape(psc(igrid)%ws(ixs^s,idir),shapes)
1022  ibuf_start=ibuf_next
1023  end do
1024  isend_r=isend_r+1
1025  call mpi_isend(sendbuffer_r(ibuf_send_r),sizes_r_send_total(i^d),&
1026  mpi_double_precision,ipe_neighbor,itag, &
1027  icomm,sendrequest_r(isend_r),ierrmpi)
1028  ibuf_send_r=ibuf_next
1029  end if
1030  end if
1031  end if
1032 
1033  end subroutine bc_send_restrict
1034 
1035  !> fill coarser neighbor's ghost cells
1036  subroutine bc_fill_restrict(igrid,i^D,iib^D)
1037  integer, intent(in) :: igrid,i^D,iib^D
1038 
1039  integer :: ic^D,n_inc^D,ixS^L,ixR^L,ipe_neighbor,ineighbor,ipole,idir
1040 
1041  ipe_neighbor=neighbor(2,i^d,igrid)
1042  if(ipe_neighbor==mype) then
1043  ic^d=1+modulo(node(pig^d_,igrid)-1,2);
1044  if({.not.(i^d==0.or.i^d==2*ic^d-3)|.or.}) return
1045  ineighbor=neighbor(1,i^d,igrid)
1046  ipole=neighbor_pole(i^d,igrid)
1047  if(ipole==0) then
1048  n_inc^d=-2*i^d+ic^d;
1049  ixs^l=ixs_r_^l(iib^d,i^d);
1050  ixr^l=ixr_r_^l(iib^d,n_inc^d);
1051  psb(ineighbor)%w(ixr^s,nwhead:nwtail)=&
1052  psc(igrid)%w(ixs^s,nwhead:nwtail)
1053  if(stagger_grid) then
1054  do idir=1,ndim
1055  ixs^l=ixs_r_stg_^l(idir,i^d);
1056  ixr^l=ixr_r_stg_^l(idir,n_inc^d);
1057  psb(ineighbor)%ws(ixr^s,idir)=psc(igrid)%ws(ixs^s,idir)
1058  end do
1059  end if
1060  else
1061  ixs^l=ixs_r_^l(iib^d,i^d);
1062  select case (ipole)
1063  {case (^d)
1064  n_inc^d=2*i^d+(3-ic^d)^d%n_inc^dd=-2*i^dd+ic^dd;\}
1065  end select
1066  ixr^l=ixr_r_^l(iib^d,n_inc^d);
1067  call pole_copy(psb(ineighbor)%w,ixg^ll,ixr^l,psc(igrid)%w,ixcog^l,ixs^l,ipole)
1068  if(stagger_grid) then
1069  do idir=1,ndim
1070  ixs^l=ixs_r_stg_^l(idir,i^d);
1071  ixr^l=ixr_r_stg_^l(idir,n_inc^d);
1072  !! Fill ghost cells
1073  call pole_copy_stg(psb(ineighbor)%ws,ixr^l,psc(igrid)%ws,ixs^l,idir,ipole)
1074  end do
1075  end if
1076  end if
1077  end if
1078 
1079  end subroutine bc_fill_restrict
1080 
1081  !> fill siblings ghost cells with received data
1082  subroutine bc_fill_srl_stg
1083  double precision :: tmp(ixGs^T)
1084  integer :: ixS^L,ixR^L,n_i^D,ixSsync^L,ixRsync^L
1085  integer :: idir
1086 
1087  ipe_neighbor=neighbor(2,i^d,igrid)
1088  if(ipe_neighbor/=mype) then
1089  ineighbor=neighbor(1,i^d,igrid)
1090  ipole=neighbor_pole(i^d,igrid)
1091 
1092  !! Now the special treatment of the pole is done here, at the receive step
1093  if (ipole==0) then
1094  ixr^l=ixr_srl_^l(iib^d,i^d);
1095  !! Unpack the buffer and fill the ghost cells
1096  n_i^d=-i^d;
1097  do idir=1,ndim
1098  ixs^l=ixs_srl_stg_^l(idir,n_i^d);
1099  ixr^l=ixr_srl_stg_^l(idir,i^d);
1100  ibuf_next=ibuf_recv_srl+sizes_srl_recv_stg(idir,i^d)
1101  tmp(ixs^s) = reshape(source=recvbuffer_srl(ibuf_recv_srl:ibuf_next-1),shape=shape(psb(igrid)%ws(ixs^s,idir)))
1102  psb(igrid)%ws(ixr^s,idir) = tmp(ixs^s)
1103  ibuf_recv_srl=ibuf_next
1104  end do
1105  else ! There is a pole
1106  select case (ipole)
1107  {case (^d)
1108  n_i^d=i^d^d%n_i^dd=-i^dd;\}
1109  end select
1110  pole_buf%ws=zero
1111  do idir=1,ndim
1112  ixr^l=ixr_srl_stg_^l(idir,i^d);
1113  ixs^l=ixs_srl_stg_^l(idir,n_i^d);
1114  ibuf_next=ibuf_recv_srl+sizes_srl_recv_stg(idir,i^d)
1115  pole_buf%ws(ixs^s,idir)=reshape(source=recvbuffer_srl(ibuf_recv_srl:ibuf_next-1),&
1116  shape=shape(psb(igrid)%ws(ixs^s,idir)))
1117  ibuf_recv_srl=ibuf_next
1118  call pole_copy_stg(psb(igrid)%ws,ixr^l,pole_buf%ws,ixs^l,idir,ipole)
1119  end do
1120  end if
1121  end if
1122 
1123  end subroutine bc_fill_srl_stg
1124 
1125  subroutine indices_for_syncing(idir,i^D,ixR^L,ixS^L,ixRsync^L,ixSsync^L)
1126  integer, intent(in) :: i^D,idir
1127  integer, intent(inout) :: ixR^L,ixS^L
1128  integer, intent(out) :: ixRsync^L,ixSsync^L
1129 
1130  ixrsync^l=ixr^l;
1131  ixssync^l=ixs^l;
1132 
1133  {
1134  if (i^d == -1 .and. idir == ^d) then
1135  ixrsyncmin^d = ixrmax^d
1136  ixrsyncmax^d = ixrmax^d
1137  ixssyncmin^d = ixsmax^d
1138  ixssyncmax^d = ixsmax^d
1139  ixrmax^d = ixrmax^d - 1
1140  ixsmax^d = ixsmax^d - 1
1141  else if (i^d == 1 .and. idir == ^d) then
1142  ixrsyncmin^d = ixrmin^d
1143  ixrsyncmax^d = ixrmin^d
1144  ixssyncmin^d = ixsmin^d
1145  ixssyncmax^d = ixsmin^d
1146  ixrmin^d = ixrmin^d + 1
1147  ixsmin^d = ixsmin^d + 1
1148  end if
1149  \}
1150 
1151  end subroutine indices_for_syncing
1152 
1153  !> fill restricted ghost cells after receipt
1155 
1156  ipole=neighbor_pole(i^d,igrid)
1157  if (ipole==0) then
1158  ! Loop over the children ic^D to and their neighbors inc^D
1159  {do ic^db=1+int((1-i^db)/2),2-int((1+i^db)/2)
1160  inc^db=2*i^db+ic^db\}
1161  ipe_neighbor=neighbor_child(2,inc^d,igrid)
1162  if(ipe_neighbor/=mype) then
1163  ineighbor=neighbor_child(1,inc^d,igrid)
1164  n_i^d=-i^d;
1165  !! Unpack the buffer and fill the ghost cells
1166  do idir=1,ndim
1167  ixr^l=ixr_r_stg_^l(idir,inc^d);
1168  ibuf_next=ibuf_recv_r+sizes_r_recv_stg(idir,inc^d)
1169  psb(igrid)%ws(ixr^s,idir)=reshape(source=recvbuffer_r(ibuf_recv_r:ibuf_next-1),&
1170  shape=shape(psb(igrid)%ws(ixr^s,idir)))
1171  ibuf_recv_r=ibuf_next
1172  end do
1173  end if
1174  {end do\}
1175  else !! There is a pole
1176  {do ic^db=1+int((1-i^db)/2),2-int((1+i^db)/2)
1177  inc^db=2*i^db+ic^db\}
1178  ipe_neighbor=neighbor_child(2,inc^d,igrid)
1179  if(ipe_neighbor/=mype) then
1180  ineighbor=neighbor_child(1,inc^d,igrid)
1181  select case(ipole)
1182  {case (^d)
1183  n_i^d=i^d^d%n_i^dd=-i^dd;\}
1184  end select
1185  ixr^l=ixr_r_^l(iib^d,inc^d);
1186  !! Unpack the buffer and fill an auxiliary array
1187  pole_buf%ws=zero
1188  do idir=1,ndim
1189  ixs^l=ixs_r_stg_^l(idir,n_i^d);
1190  ixr^l=ixr_r_stg_^l(idir,inc^d);
1191  ibuf_next=ibuf_recv_r+sizes_r_recv_stg(idir,inc^d)
1192  pole_buf%ws(ixr^s,idir)=reshape(source=recvbuffer_r(ibuf_recv_r:ibuf_next-1),&
1193  shape=shape(psb(igrid)%ws(ixr^s,idir)))
1194  call pole_copy_stg(psb(igrid)%ws,ixr^l,pole_buf%ws,ixr^l,idir,ipole)
1195  ibuf_recv_r=ibuf_next
1196  end do
1197  end if
1198  {end do\}
1199  end if
1200 
1201  end subroutine bc_fill_restrict_stg
1202 
1203  !> Receive from coarse neighbor
1204  subroutine bc_recv_prolong
1205 
1206  ic^d=1+modulo(node(pig^d_,igrid)-1,2);
1207  if ({.not.(i^d==0.or.i^d==2*ic^d-3)|.or.}) return
1208 
1209  ipe_neighbor=neighbor(2,i^d,igrid)
1210  if (ipe_neighbor/=mype) then
1211  irecv_c=irecv_c+1
1212  inc^d=ic^d+i^d;
1213  itag=(3**^nd+4**^nd)*(igrid-1)+3**^nd+{inc^d*4**(^d-1)+}
1214  call mpi_irecv(psc(igrid)%w,1,type_recv_p(iib^d,inc^d), &
1215  ipe_neighbor,itag,icomm,recvrequest_c_p(irecv_c),ierrmpi)
1216  if(stagger_grid) then
1217  irecv_p=irecv_p+1
1218  call mpi_irecv(recvbuffer_p(ibuf_recv_p),sizes_p_recv_total(inc^d),&
1219  mpi_double_precision,ipe_neighbor,itag,&
1220  icomm,recvrequest_p(irecv_p),ierrmpi)
1222  end if
1223  end if
1224 
1225  end subroutine bc_recv_prolong
1226 
1227  !> Send to finer neighbor
1228  subroutine bc_send_prolong
1229  integer :: ii^D
1230 
1231  ipole=neighbor_pole(i^d,igrid)
1232 
1233  {do ic^db=1+int((1-i^db)/2),2-int((1+i^db)/2)
1234  inc^db=2*i^db+ic^db\}
1235  ipe_neighbor=neighbor_child(2,inc^d,igrid)
1236  if(ipe_neighbor/=mype) then
1237  ixs^l=ixs_p_^l(iib^d,inc^d);
1238  ineighbor=neighbor_child(1,inc^d,igrid)
1239  if(ipole==0) then
1240  n_i^d=-i^d;
1241  n_inc^d=ic^d+n_i^d;
1242  isend_c=isend_c+1
1243  itag=(3**^nd+4**^nd)*(ineighbor-1)+3**^nd+{n_inc^d*4**(^d-1)+}
1244  call mpi_isend(psb(igrid)%w,1,type_send_p(iib^d,inc^d), &
1245  ipe_neighbor,itag,icomm,sendrequest_c_p(isend_c),ierrmpi)
1246  if(stagger_grid) then
1247  ibuf_start=ibuf_send_p
1248  do idir=1,ndim
1249  ixs^l=ixs_p_stg_^l(idir,inc^d);
1250  ibuf_next=ibuf_start+sizes_p_send_stg(idir,inc^d)
1251  shapes=(/sizes_p_send_stg(idir,inc^d)/)
1252  sendbuffer_p(ibuf_start:ibuf_next-1)=&
1253  reshape(psb(igrid)%ws(ixs^s,idir),shapes)
1254  ibuf_start=ibuf_next
1255  end do
1256  isend_p=isend_p+1
1257  call mpi_isend(sendbuffer_p(ibuf_send_p),sizes_p_send_total(inc^d),&
1258  mpi_double_precision,ipe_neighbor,itag, &
1259  icomm,sendrequest_p(isend_p),ierrmpi)
1260  ibuf_send_p=ibuf_next
1261  end if
1262  else
1263  select case (ipole)
1264  {case (^d)
1265  n_inc^d=inc^d^d%n_inc^dd=ic^dd-i^dd;\}
1266  end select
1267  if(isend_buf(ipwbuf)/=0) then
1268  call mpi_wait(sendrequest_c_p(isend_buf(ipwbuf)), &
1269  sendstatus_c_p(:,isend_buf(ipwbuf)),ierrmpi)
1270  deallocate(pwbuf(ipwbuf)%w)
1271  end if
1272  allocate(pwbuf(ipwbuf)%w(ixs^s,nwhead:nwtail))
1273  call pole_buffer(pwbuf(ipwbuf)%w,ixs^l,ixs^l,psb(igrid)%w,ixg^ll,ixs^l)
1274  isend_c=isend_c+1
1275  isend_buf(ipwbuf)=isend_c
1276  itag=(3**^nd+4**^nd)*(ineighbor-1)+3**^nd+{n_inc^d*4**(^d-1)+}
1277  isizes={(ixsmax^d-ixsmin^d+1)*}*nwbc
1278  call mpi_isend(pwbuf(ipwbuf)%w,isizes,mpi_double_precision, &
1279  ipe_neighbor,itag,icomm,sendrequest_c_p(isend_c),ierrmpi)
1280  ipwbuf=1+modulo(ipwbuf,npwbuf)
1281  if(stagger_grid) then
1282  ibuf_start=ibuf_send_p
1283  do idir=1,ndim
1284  ixs^l=ixs_p_stg_^l(idir,inc^d);
1285  ibuf_next=ibuf_start+sizes_p_send_stg(idir,inc^d)
1286  shapes=(/sizes_p_send_stg(idir,inc^d)/)
1287  sendbuffer_p(ibuf_start:ibuf_next-1)=&
1288  reshape(psb(igrid)%ws(ixs^s,idir),shapes)
1289  ibuf_start=ibuf_next
1290  end do
1291  isend_p=isend_p+1
1292  call mpi_isend(sendbuffer_p(ibuf_send_p),sizes_p_send_total(inc^d),&
1293  mpi_double_precision,ipe_neighbor,itag, &
1294  icomm,sendrequest_p(isend_p),ierrmpi)
1295  ibuf_send_p=ibuf_next
1296  end if
1297  end if
1298  end if
1299  {end do\}
1300 
1301  end subroutine bc_send_prolong
1302 
1303  !> Send to finer neighbor
1304  subroutine bc_fill_prolong(igrid,i^D,iib^D)
1305  integer, intent(in) :: igrid,i^D,iib^D
1306 
1307  integer :: ipe_neighbor,ineighbor,ixS^L,ixR^L,ic^D,inc^D,ipole,idir
1308 
1309  ipole=neighbor_pole(i^d,igrid)
1310 
1311  if(ipole==0) then
1312  {do ic^db=1+int((1-i^db)/2),2-int((1+i^db)/2)
1313  inc^db=2*i^db+ic^db\}
1314  ipe_neighbor=neighbor_child(2,inc^d,igrid)
1315  if(ipe_neighbor==mype) then
1316  ixs^l=ixs_p_^l(iib^d,inc^d);
1317  ineighbor=neighbor_child(1,inc^d,igrid)
1318  ipole=neighbor_pole(i^d,igrid)
1319  n_i^d=-i^d;
1320  n_inc^d=ic^d+n_i^d;
1321  ixr^l=ixr_p_^l(iib^d,n_inc^d);
1322  psc(ineighbor)%w(ixr^s,nwhead:nwtail) &
1323  =psb(igrid)%w(ixs^s,nwhead:nwtail)
1324  if(stagger_grid) then
1325  do idir=1,ndim
1326  ixs^l=ixs_p_stg_^l(idir,inc^d);
1327  ixr^l=ixr_p_stg_^l(idir,n_inc^d);
1328  psc(ineighbor)%ws(ixr^s,idir)=psb(igrid)%ws(ixs^s,idir)
1329  end do
1330  end if
1331  end if
1332  {end do\}
1333  else
1334  {do ic^db=1+int((1-i^db)/2),2-int((1+i^db)/2)
1335  inc^db=2*i^db+ic^db\}
1336  ipe_neighbor=neighbor_child(2,inc^d,igrid)
1337  if(ipe_neighbor==mype) then
1338  ixs^l=ixs_p_^l(iib^d,inc^d);
1339  ineighbor=neighbor_child(1,inc^d,igrid)
1340  ipole=neighbor_pole(i^d,igrid)
1341  select case (ipole)
1342  {case (^d)
1343  n_inc^d=inc^d^d%n_inc^dd=ic^dd-i^dd;\}
1344  end select
1345  ixr^l=ixr_p_^l(iib^d,n_inc^d);
1346  call pole_copy(psc(ineighbor)%w,ixcog^l,ixr^l,psb(igrid)%w,ixg^ll,ixs^l,ipole)
1347  if(stagger_grid) then
1348  do idir=1,ndim
1349  ixs^l=ixs_p_stg_^l(idir,inc^d);
1350  ixr^l=ixr_p_stg_^l(idir,n_inc^d);
1351  call pole_copy_stg(psc(ineighbor)%ws,ixr^l,psb(igrid)%ws,ixs^l,idir,ipole)
1352  end do
1353  end if
1354  end if
1355  {end do\}
1356  end if
1357  end subroutine bc_fill_prolong
1358 
1359  subroutine gc_prolong(igrid)
1360  integer, intent(in) :: igrid
1361 
1362  integer :: iib^D,i^D,idims,iside
1363  logical,dimension(-1:1^D&) :: NeedProlong
1364 
1365  ^d&iib^d=idphyb(^d,igrid);
1366  needprolong=.false.
1367  {do i^db=-1,1\}
1368  if (skip_direction([ i^d ])) cycle
1369  if (neighbor_type(i^d,igrid)==neighbor_coarse) then
1370  call bc_prolong(igrid,i^d,iib^d)
1371  needprolong(i^d)=.true.
1372  end if
1373  {end do\}
1374  if(stagger_grid) then
1375  ! Ghost cell prolongation for staggered variables
1376  ! must be done in a specific order.
1377  ! First the first neighbours, which have 2 indices=0 in 3D
1378  ! or one index=0 in 2D
1379  block=>psb(igrid)
1380  do idims=1,ndim
1381  i^d=0;
1382  select case(idims)
1383  {case(^d)
1384  do i^d=-1,1,2
1385  if (needprolong(i^dd)) call bc_prolong_stg(igrid,i^dd,iib^dd,needprolong)
1386  end do
1387  \}
1388  end select
1389  end do
1390  ! Then the second neighbours which have 1 index=0 in 3D
1391  ! (Only in 3D)
1392  {^ifthreed
1393  i1=0;
1394  do i2=-1,1,2
1395  do i3=-1,1,2
1396  if (needprolong(i^d)) call bc_prolong_stg(igrid,i^d,iib^d,needprolong)
1397  end do
1398  end do
1399  i2=0;
1400  do i3=-1,1,2
1401  do i1=-1,1,2
1402  if (needprolong(i^d)) call bc_prolong_stg(igrid,i^d,iib^d,needprolong)
1403  end do
1404  end do
1405  i3=0;
1406  do i1=-1,1,2
1407  do i2=-1,1,2
1408  if (needprolong(i^d)) call bc_prolong_stg(igrid,i^d,iib^d,needprolong)
1409  end do
1410  end do
1411  }
1412  ! Finally, the corners, that have no index=0
1413  {do i^d=-1,1,2\}
1414  if (needprolong(i^d)) call bc_prolong_stg(igrid,i^d,iib^d,needprolong)
1415  {end do\}
1416  end if
1417  end subroutine gc_prolong
1418 
1419  !> fill coarser representative with data from coarser neighbors
1421  ic^d=1+modulo(node(pig^d_,igrid)-1,2);
1422  if ({.not.(i^d==0.or.i^d==2*ic^d-3)|.or.}) return
1423 
1424  ipe_neighbor=neighbor(2,i^d,igrid)
1425  if(ipe_neighbor/=mype) then
1426  ineighbor=neighbor(1,i^d,igrid)
1427  ipole=neighbor_pole(i^d,igrid)
1428 
1429  if (ipole==0) then !! There is no pole
1430  inc^d=ic^d+i^d;
1431  ixr^l=ixr_p_^l(iib^d,inc^d);
1432  do idir=1,ndim
1433  ixr^l=ixr_p_stg_^l(idir,inc^d);
1434  ibuf_next=ibuf_recv_p+sizes_p_recv_stg(idir,inc^d)
1435  psc(igrid)%ws(ixr^s,idir)=reshape(source=recvbuffer_p(ibuf_recv_p:ibuf_next-1),&
1436  shape=shape(psc(igrid)%ws(ixr^s,idir)))
1437  ibuf_recv_p=ibuf_next
1438  end do
1439  else !! There is a pole
1440  inc^d=ic^d+i^d;
1441  select case (ipole)
1442  {case (^d)
1443  n_inc^d=2*i^d+(3-ic^d)^d%n_inc^dd=-2*i^dd+ic^dd;\}
1444  end select
1445  !! Unpack the buffer and fill an auxiliary array
1446  pole_buf%ws=zero
1447  do idir=1,ndim
1448  ixr^l=ixr_p_stg_^l(idir,inc^d);
1449  ibuf_next=ibuf_recv_p+sizes_p_recv_stg(idir,inc^d)
1450  pole_buf%ws(ixr^s,idir)=reshape(source=recvbuffer_p(ibuf_recv_p:ibuf_next-1),&
1451  shape=shape(psc(igrid)%ws(ixr^s,idir)))
1452  call pole_copy_stg(psc(igrid)%ws,ixr^l,pole_buf%ws,ixr^l,idir,ipole)
1453  ibuf_recv_p=ibuf_next
1454  end do
1455  end if
1456  end if
1457 
1458  end subroutine bc_fill_prolong_stg
1459 
1460  !> do prolongation for fine blocks after receipt data from coarse neighbors
1461  subroutine bc_prolong(igrid,i^D,iib^D)
1463 
1464  integer :: i^D,iib^D,igrid
1465  integer :: ixFi^L,ixCo^L,ii^D, idims,iside,ixB^L
1466  double precision :: dxFi^D, dxCo^D, xFimin^D, xComin^D, invdxCo^D
1467 
1468  ixfi^l=ixr_srl_^l(iib^d,i^d);
1469  dxfi^d=rnode(rpdx^d_,igrid);
1470  dxco^d=two*dxfi^d;
1471  invdxco^d=1.d0/dxco^d;
1472 
1473  ! compute the enlarged grid lower left corner coordinates
1474  ! these are true coordinates for an equidistant grid,
1475  ! but we can temporarily also use them for getting indices
1476  ! in stretched grids
1477  xfimin^d=rnode(rpxmin^d_,igrid)-dble(nghostcells)*dxfi^d;
1478  xcomin^d=rnode(rpxmin^d_,igrid)-dble(nghostcells)*dxco^d;
1479 
1480  if(stagger_grid.and.phyboundblock(igrid).and.bcphys) then
1481  block=>psc(igrid)
1482  do idims=1,ndim
1483  ixcomin^d=int((xfimin^d+(dble(ixfimin^d)-half)*dxfi^d-xcomin^d)*invdxco^d)+1-1;
1484  ixcomax^d=int((xfimin^d+(dble(ixfimax^d)-half)*dxfi^d-xcomin^d)*invdxco^d)+1+1;
1485  {^ifthreed
1486  ! avoid using undetermined ghost cells at physical boundary edges
1487  if(idims == 1) then
1488  if(neighbor_type(-1,0,0,igrid)==neighbor_boundary .or. &
1489  neighbor_type(1,0,0,igrid)==neighbor_boundary) then
1490  if(neighbor_type(0,-1,0,igrid)==neighbor_boundary) ixcomin2=ixcommin2
1491  if(neighbor_type(0,0,-1,igrid)==neighbor_boundary) ixcomin3=ixcommin3
1492  if(neighbor_type(0,1,0,igrid)==neighbor_boundary) ixcomax2=ixcommax2
1493  if(neighbor_type(0,0,1,igrid)==neighbor_boundary) ixcomax3=ixcommax3
1494  end if
1495  else if(idims == 2) then
1496  if(neighbor_type(0,-1,0,igrid)==neighbor_boundary .or. &
1497  neighbor_type(0,1,0,igrid)==neighbor_boundary) then
1498  if(neighbor_type(0,0,-1,igrid)==neighbor_boundary) ixcomin3=ixcommin3
1499  if(neighbor_type(0,0,1,igrid)==neighbor_boundary) ixcomax3=ixcommax3
1500  end if
1501  end if
1502  }
1503  do iside=1,2
1504  ii^d=kr(^d,idims)*(2*iside-3);
1505  if(neighbor_type(ii^d,igrid)/=neighbor_boundary) cycle
1506  if(( {(iside==1.and.idims==^d.and.ixcomin^d<ixcogmin^d+nghostcells)|.or.} ) &
1507  .or.( {(iside==2.and.idims==^d.and.ixcomax^d>ixcogmax^d-nghostcells)|.or. })) then
1508  {ixbmin^d=merge(ixcogmin^d,ixcomin^d,idims==^d);}
1509  {ixbmax^d=merge(ixcogmax^d,ixcomax^d,idims==^d);}
1510  call bc_phys(iside,idims,time,0.d0,psc(igrid),ixcog^l,ixb^l)
1511  end if
1512  end do
1513  end do
1514  end if
1515 
1516  if(prolongprimitive) then
1517  ! following line again assumes equidistant grid, but
1518  ! just computes indices, so also ok for stretched case
1519  ! reason for +1-1 and +1+1: the coarse representation has
1520  ! also nghostcells at each side. During
1521  ! prolongation, we need cells to left and right, hence -1/+1
1522  block=>psc(igrid)
1523  ixcomin^d=int((xfimin^d+(dble(ixfimin^d)-half)*dxfi^d-xcomin^d)*invdxco^d)+1-1;
1524  ixcomax^d=int((xfimin^d+(dble(ixfimax^d)-half)*dxfi^d-xcomin^d)*invdxco^d)+1+1;
1525  call phys_to_primitive(ixcog^l,ixco^l,&
1526  psc(igrid)%w,psc(igrid)%x)
1527  end if
1528 
1529  if(ghost_copy) then
1530  call interpolation_copy(igrid,ixfi^l,dxfi^d,xfimin^d,dxco^d,invdxco^d,xcomin^d)
1531  else
1532  call interpolation_linear(igrid,ixfi^l,dxfi^d,xfimin^d,dxco^d,invdxco^d,xcomin^d)
1533  end if
1534 
1535  if(prolongprimitive) then
1536  block=>psc(igrid)
1537  call phys_to_conserved(ixcog^l,ixco^l,psc(igrid)%w,psc(igrid)%x)
1538  end if
1539 
1540  end subroutine bc_prolong
1541 
1542  subroutine bc_prolong_stg(igrid,i^D,iib^D,NeedProlong)
1543  use mod_amr_fct
1544  integer :: igrid,i^D,iib^D
1545  logical,dimension(-1:1^D&) :: NeedProlong
1546  logical :: fine_^Lin
1547  integer :: ixFi^L,ixCo^L
1548  double precision :: dxFi^D,dxCo^D,xFimin^D,xComin^D,invdxCo^D
1549  ! Check what is already at the desired level
1550  fine_^lin=.false.;
1551  {
1552  if(i^d>-1) fine_min^din=(.not.needprolong(i^dd-kr(^d,^dd)).and.neighbor_type(i^dd-kr(^d,^dd),igrid)/=1)
1553  if(i^d<1) fine_max^din=(.not.needprolong(i^dd+kr(^d,^dd)).and.neighbor_type(i^dd+kr(^d,^dd),igrid)/=1)
1554  \}
1555 
1556  ixfi^l=ixr_srl_^l(iib^d,i^d);
1557 
1558  dxfi^d=rnode(rpdx^d_,igrid);
1559  dxco^d=two*dxfi^d;
1560  invdxco^d=1.d0/dxco^d;
1561 
1562  xfimin^d=rnode(rpxmin^d_,igrid)-dble(nghostcells)*dxfi^d;
1563  xcomin^d=rnode(rpxmin^d_,igrid)-dble(nghostcells)*dxco^d;
1564 
1565  ! moved the physical boundary filling here, to only fill the
1566  ! part needed
1567 
1568  ixcomin^d=int((xfimin^d+(dble(ixfimin^d)-half)*dxfi^d-xcomin^d)*invdxco^d)+1-1;
1569  ixcomax^d=int((xfimin^d+(dble(ixfimax^d)-half)*dxfi^d-xcomin^d)*invdxco^d)+1+1;
1570 
1571  if(prolongprimitive) call phys_to_primitive(ixg^ll,ixfi^l,psb(igrid)%w,psb(igrid)%x)
1572 
1573  call prolong_2nd_stg(psc(igrid),psb(igrid),ixco^l,ixfi^l,dxco^d,xcomin^d,dxfi^d,xfimin^d,.true.,fine_^lin)
1574 
1575  if(prolongprimitive) call phys_to_conserved(ixg^ll,ixfi^l,psb(igrid)%w,psb(igrid)%x)
1576 
1577  ! The current region has already been refined, so it does not need to be prolonged again
1578  needprolong(i^d)=.false.
1579 
1580  end subroutine bc_prolong_stg
1581 
1582  subroutine interpolation_linear(igrid,ixFi^L,dxFi^D,xFimin^D, &
1583  dxCo^D,invdxCo^D,xComin^D)
1584  use mod_physics, only: phys_to_conserved
1585  integer, intent(in) :: igrid, ixFi^L
1586  double precision, intent(in) :: dxFi^D, xFimin^D,dxCo^D, invdxCo^D, xComin^D
1587 
1588  integer :: ixCo^D, jxCo^D, hxCo^D, ixFi^D, ix^D, iw, idims, nwmin,nwmax
1589  double precision :: xCo^D, xFi^D, eta^D
1590  double precision :: slopeL, slopeR, slopeC, signC, signR
1591  double precision :: slope(1:nw,ndim)
1592  !!double precision :: local_invdxCo^D
1593  double precision :: signedfactorhalf^D
1594  !integer :: ixshift^D, icase
1595 
1596  !icase=mod(nghostcells,2)
1597 
1598  if(prolongprimitive) then
1599  nwmin=1
1600  nwmax=nw
1601  else
1602  nwmin=nwhead
1603  nwmax=nwtail
1604  end if
1605 
1606  {do ixfi^db = ixfi^lim^db
1607  ! cell-centered coordinates of fine grid point
1608  ! here we temporarily use an equidistant grid
1609  xfi^db=xfimin^db+(dble(ixfi^db)-half)*dxfi^db
1610 
1611  ! indices of coarse cell which contains the fine cell
1612  ! since we computed lower left corner earlier
1613  ! in equidistant fashion: also ok for stretched case
1614  ixco^db=int((xfi^db-xcomin^db)*invdxco^db)+1
1615 
1616  ! cell-centered coordinates of coarse grid point
1617  ! here we temporarily use an equidistant grid
1618  xco^db=xcomin^db+(dble(ixco^db)-half)*dxco^db \}
1619 
1620  !if(.not.slab) then
1621  ! ^D&local_invdxCo^D=1.d0/psc(igrid)%dx({ixCo^DD},^D)\
1622  !endif
1623 
1624  if(slab_uniform) then
1625  ! actual cell-centered coordinates of fine grid point
1626  !!^D&xFi^D=block%x({ixFi^DD},^D)\
1627  ! actual cell-centered coordinates of coarse grid point
1628  !!^D&xCo^D=psc(igrid)%x({ixCo^DD},^D)\
1629  ! normalized distance between fine/coarse cell center
1630  ! in coarse cell: ranges from -0.5 to 0.5 in each direction
1631  ! (origin is coarse cell center)
1632  ! this is essentially +1/4 or -1/4 on cartesian mesh
1633  eta^d=(xfi^d-xco^d)*invdxco^d;
1634  else
1635  !select case(icase)
1636  ! case(0)
1637  !{! here we assume an even number of ghostcells!!!
1638  !ixshift^D=2*(mod(ixFi^D,2)-1)+1
1639  !if(ixshift^D>0.0d0)then
1640  ! ! oneven fine grid points
1641  ! eta^D=-0.5d0*(one-block%dvolume(ixFi^DD) &
1642  ! /sum(block%dvolume(ixFi^D:ixFi^D+1^D%ixFi^DD)))
1643  !else
1644  ! ! even fine grid points
1645  ! eta^D=+0.5d0*(one-block%dvolume(ixFi^DD) &
1646  ! /sum(block%dvolume(ixFi^D-1:ixFi^D^D%ixFi^DD)))
1647  !endif\}
1648  ! case(1)
1649  !{! here we assume an odd number of ghostcells!!!
1650  !ixshift^D=2*(mod(ixFi^D,2)-1)+1
1651  !if(ixshift^D>0.0d0)then
1652  ! ! oneven fine grid points
1653  ! eta^D=+0.5d0*(one-block%dvolume(ixFi^DD) &
1654  ! /sum(block%dvolume(ixFi^D-1:ixFi^D^D%ixFi^DD)))
1655  !else
1656  ! ! even fine grid points
1657  ! eta^D=-0.5d0*(one-block%dvolume(ixFi^DD) &
1658  ! /sum(block%dvolume(ixFi^D:ixFi^D+1^D%ixFi^DD)))
1659  !endif\}
1660  ! case default
1661  ! call mpistop("no such case")
1662  !end select
1663  ! the different cases for even/uneven number of ghost cells
1664  ! are automatically handled using the relative index to ixMlo
1665  ! as well as the pseudo-coordinates xFi and xCo
1666  ! these latter differ from actual cell centers when stretching is used
1667  ix^d=2*int((ixfi^d+ixmlo^d)/2)-ixmlo^d;
1668  {if(xfi^d>xco^d) then
1669  signedfactorhalf^d=0.5d0
1670  else
1671  signedfactorhalf^d=-0.5d0
1672  end if
1673  eta^d=signedfactorhalf^d*(one-psb(igrid)%dvolume(ixfi^dd) &
1674  /sum(psb(igrid)%dvolume(ix^d:ix^d+1^d%ixFi^dd))) \}
1675  !{eta^D=(xFi^D-xCo^D)*invdxCo^D &
1676  ! *two*(one-block%dvolume(ixFi^DD) &
1677  ! /sum(block%dvolume(ix^D:ix^D+1^D%ixFi^DD))) \}
1678  end if
1679 
1680  do idims=1,ndim
1681  hxco^d=ixco^d-kr(^d,idims)\
1682  jxco^d=ixco^d+kr(^d,idims)\
1683 
1684  do iw=nwmin,nwmax
1685  slopel=psc(igrid)%w(ixco^d,iw)-psc(igrid)%w(hxco^d,iw)
1686  sloper=psc(igrid)%w(jxco^d,iw)-psc(igrid)%w(ixco^d,iw)
1687  slopec=half*(sloper+slopel)
1688 
1689  ! get limited slope
1690  signr=sign(one,sloper)
1691  signc=sign(one,slopec)
1692  !select case(prolong_limiter)
1693  !case(1)
1694  ! ! unlimit
1695  ! slope(iw,idims)=slopeC
1696  !case(2)
1697  ! ! minmod
1698  ! slope(iw,idims)=signR*max(zero,min(dabs(slopeR), &
1699  ! signR*slopeL))
1700  !case(3)
1701  ! ! woodward
1702  ! slope(iw,idims)=two*signR*max(zero,min(dabs(slopeR), &
1703  ! signR*slopeL,signR*half*slopeC))
1704  !case(4)
1705  ! ! koren
1706  ! slope(iw,idims)=signR*max(zero,min(two*signR*slopeL, &
1707  ! (dabs(slopeR)+two*slopeL*signR)*third,two*dabs(slopeR)))
1708  !case default
1709  slope(iw,idims)=signc*max(zero,min(dabs(slopec), &
1710  signc*slopel,signc*sloper))
1711  !end select
1712  end do
1713  end do
1714 
1715  ! Interpolate from coarse cell using limited slopes
1716  psb(igrid)%w(ixfi^d,nwmin:nwmax)=psc(igrid)%w(ixco^d,nwmin:nwmax)+&
1717  {(slope(nwmin:nwmax,^d)*eta^d)+}
1718 
1719  {end do\}
1720 
1721  if(prolongprimitive) then
1722  block=>psb(igrid)
1723  call phys_to_conserved(ixg^ll,ixfi^l,psb(igrid)%w,psb(igrid)%x)
1724  end if
1725 
1726  end subroutine interpolation_linear
1727 
1728  subroutine interpolation_copy(igrid, ixFi^L,dxFi^D,xFimin^D, &
1729  dxCo^D,invdxCo^D,xComin^D)
1730  use mod_physics, only: phys_to_conserved
1731  integer, intent(in) :: igrid, ixFi^L
1732  double precision, intent(in) :: dxFi^D, xFimin^D,dxCo^D, invdxCo^D, xComin^D
1733 
1734  integer :: ixCo^D, ixFi^D, nwmin,nwmax
1735  double precision :: xFi^D
1736 
1737  if(prolongprimitive) then
1738  nwmin=1
1739  nwmax=nw
1740  else
1741  nwmin=nwhead
1742  nwmax=nwtail
1743  end if
1744 
1745  {do ixfi^db = ixfi^lim^db
1746  ! cell-centered coordinates of fine grid point
1747  xfi^db=xfimin^db+(dble(ixfi^db)-half)*dxfi^db
1748 
1749  ! indices of coarse cell which contains the fine cell
1750  ! note: this also works for stretched grids
1751  ixco^db=int((xfi^db-xcomin^db)*invdxco^db)+1\}
1752 
1753  ! Copy from coarse cell
1754  psb(igrid)%w(ixfi^d,nwmin:nwmax)=psc(igrid)%w(ixco^d,nwmin:nwmax)
1755 
1756  {end do\}
1757 
1758  if(prolongprimitive) call phys_to_conserved(ixg^ll,ixfi^l,psb(igrid)%w,psb(igrid)%x)
1759 
1760  end subroutine interpolation_copy
1761 
1762  subroutine pole_copy(wrecv,ixIR^L,ixR^L,wsend,ixIS^L,ixS^L,ipole)
1763 
1764  integer, intent(in) :: ixIR^L,ixR^L,ixIS^L,ixS^L,ipole
1765  double precision :: wrecv(ixIR^S,1:nw), wsend(ixIS^S,1:nw)
1766 
1767  integer :: iw, iside, iB
1768 
1769  select case (ipole)
1770  {case (^d)
1771  iside=int((i^d+3)/2)
1772  ib=2*(^d-1)+iside
1773  do iw=nwhead,nwtail
1774  select case (typeboundary(iw,ib))
1775  case (bc_symm)
1776  wrecv(ixr^s,iw) = wsend(ixsmax^d:ixsmin^d:-1^d%ixS^s,iw)
1777  case (bc_asymm)
1778  wrecv(ixr^s,iw) =-wsend(ixsmax^d:ixsmin^d:-1^d%ixS^s,iw)
1779  case default
1780  call mpistop("Pole boundary condition should be symm or asymm")
1781  end select
1782  end do \}
1783  end select
1784 
1785  end subroutine pole_copy
1786 
1787  subroutine pole_copy_stg(wrecv,ixR^L,wsend,ixS^L,idirs,ipole)
1788 
1789  integer, intent(in) :: ixR^L,ixS^L,idirs,ipole
1790 
1791  double precision :: wrecv(ixGs^T,1:nws), wsend(ixGs^T,1:nws)
1792  integer :: iB, iside
1793 
1794  select case (ipole)
1795  {case (^d)
1796  iside=int((i^d+3)/2)
1797  ib=2*(^d-1)+iside
1798  select case (typeboundary(iw_mag(idirs),ib))
1799  case (bc_symm)
1800  wrecv(ixr^s,idirs) = wsend(ixsmax^d:ixsmin^d:-1^d%ixS^s,idirs)
1801  case (bc_asymm)
1802  wrecv(ixr^s,idirs) =-wsend(ixsmax^d:ixsmin^d:-1^d%ixS^s,idirs)
1803  case default
1804  call mpistop("Pole boundary condition should be symm or asymm")
1805  end select
1806  \}
1807  end select
1808 
1809  end subroutine pole_copy_stg
1810 
1811  subroutine pole_buffer(wrecv,ixIR^L,ixR^L,wsend,ixIS^L,ixS^L)
1812 
1813  integer, intent(in) :: ixIR^L,ixR^L,ixIS^L,ixS^L
1814  double precision :: wrecv(ixIR^S,nwhead:nwtail), wsend(ixIS^S,1:nw)
1815 
1816  integer :: iw, iside, iB
1817 
1818  select case (ipole)
1819  {case (^d)
1820  iside=int((i^d+3)/2)
1821  ib=2*(^d-1)+iside
1822  do iw=nwhead,nwtail
1823  select case (typeboundary(iw,ib))
1824  case (bc_symm)
1825  wrecv(ixr^s,iw) = wsend(ixsmax^d:ixsmin^d:-1^d%ixS^s,iw)
1826  case (bc_asymm)
1827  wrecv(ixr^s,iw) =-wsend(ixsmax^d:ixsmin^d:-1^d%ixS^s,iw)
1828  case default
1829  call mpistop("Pole boundary condition should be symm or asymm")
1830  end select
1831  end do \}
1832  end select
1833 
1834  end subroutine pole_buffer
1835 
1836  end subroutine getbc
1837 
1838  subroutine identifyphysbound(s,iib^D)
1840 
1841  type(state) :: s
1842  integer, intent(out) :: iib^D
1843 
1844  {
1845  if(s%is_physical_boundary(2*^d) .and. &
1846  s%is_physical_boundary(2*^d-1)) then
1847  iib^d=2
1848  else if(s%is_physical_boundary(2*^d-1)) then
1849  iib^d=-1
1850  else if(s%is_physical_boundary(2*^d)) then
1851  iib^d=1
1852  else
1853  iib^d=0
1854  end if
1855  \}
1856 
1857  end subroutine identifyphysbound
1858 
1859 end module mod_ghostcells_update
subroutine bc_phys(iside, idims, time, qdt, s, ixGL, ixBL)
fill ghost cells at a physical boundary
subroutine getintbc(time, ixGL)
fill inner boundary values
subroutine coarsen_grid(sFi, ixFiGL, ixFiL, sCo, ixCoGL, ixCoL)
coarsen one grid to its coarser representative
Definition: coarsen.t:3
subroutine mpistop(message)
Exit MPI-AMRVAC with an error message.
Definition: comm_lib.t:194
subroutine bc_recv_restrict
Receive from fine neighbor.
subroutine gc_prolong(igrid)
subroutine interpolation_linear(igrid, ixFiL, dxFiD, xFiminD, dxCoD, invdxCoD, xCominD)
subroutine bc_prolong(igrid, iD, iibD)
do prolongation for fine blocks after receipt data from coarse neighbors
subroutine indices_for_syncing(idir, iD, ixRL, ixSL, ixRsyncL, ixSsyncL)
subroutine bc_fill_srl_stg
fill siblings ghost cells with received data
subroutine pole_copy(wrecv, ixIRL, ixRL, wsend, ixISL, ixSL, ipole)
subroutine pole_copy_stg(wrecv, ixRL, wsend, ixSL, idirs, ipole)
subroutine bc_recv_prolong
Receive from coarse neighbor.
subroutine bc_fill_prolong(igrid, iD, iibD)
Send to finer neighbor.
subroutine bc_recv_srl
Receive from sibling at same refinement level.
subroutine bc_fill_restrict(igrid, iD, iibD)
fill coarser neighbor's ghost cells
subroutine fill_coarse_boundary(igrid, iD)
subroutine interpolation_copy(igrid, ixFiL, dxFiD, xFiminD, dxCoD, invdxCoD, xCominD)
subroutine bc_send_restrict
Send to coarser neighbor.
subroutine pole_buffer(wrecv, ixIRL, ixRL, wsend, ixISL, ixSL)
subroutine fill_boundary_before_gc(igrid)
Physical boundary conditions.
subroutine bc_fill_srl(igrid, iD, iibD)
subroutine bc_fill_restrict_stg
fill restricted ghost cells after receipt
subroutine bc_send_srl
Send to sibling at same refinement level.
subroutine fill_boundary_after_gc(igrid)
Physical boundary conditions.
logical function skip_direction(dir)
subroutine bc_send_prolong
Send to finer neighbor.
subroutine bc_fill_prolong_stg
fill coarser representative with data from coarser neighbors
subroutine bc_prolong_stg(igrid, iD, iibD, NeedProlong)
subroutine, public prolong_2nd_stg(sCo, sFi, ixCoLin, ixFiLin, dxCoD, xCominD, dxFiD, xFiminD, ghost, fine_Lin)
This subroutine performs a 2nd order prolongation for a staggered field F, preserving the divergence ...
Definition: mod_amr_fct.t:41
update ghost cells of all blocks including physical boundaries
integer, dimension(-1:2^d &,-1:1^d &), target type_recv_srl_f
integer, dimension(-1:1^d &) sizes_r_send_total
integer, dimension(-1:1^d &, 0:3^d &), target type_recv_p_p1
integer, dimension(-1:2^d &,-1:1^d &), target type_recv_srl_p2
integer, dimension(-1:1^d &,-1:1^d &), target type_send_r_p2
integer, dimension(^nd, 0:3) l
integer, dimension(:), allocatable sendrequest_c_p
integer, dimension(^nd,-1:1) ixs_r_stg_
integer, dimension(^nd,-1:1) ixs_srl_stg_
subroutine get_bc_comm_type(comm_type, ixL, ixGL, nwstart, nwbc)
integer, dimension(^nd, 0:3^d &) sizes_p_send_stg
integer, dimension(-1:1^d &) sizes_srl_send_total
integer, dimension(-1:1^d &, 0:3^d &), target type_send_p_p1
subroutine identifyphysbound(s, iibD)
integer, dimension(^nd, 0:3) ixr_r_stg_
integer, dimension(:), allocatable sendrequest_r
integer, dimension(-1:2^d &,-1:1^d &), target type_send_srl_p1
integer, dimension(:), allocatable recvrequest_r
double precision, dimension(:), allocatable sendbuffer_p
integer, dimension(:,:), allocatable sendstatus_c_p
integer, dimension(^nd,-1:1^d &) sizes_r_send_stg
integer, dimension(-1:1^d &, 0:3^d &), target type_send_p_f
integer, dimension(-1:2,-1:1) ixr_srl_
integer, dimension(:^d &,:^d &), pointer type_recv_srl
integer, dimension(:,:), allocatable sendstatus_c_sr
integer, dimension(-1:1^d &, 0:3^d &), target type_recv_r_p2
double precision, dimension(:), allocatable recvbuffer_r
integer, dimension(:), allocatable recvrequest_c_p
subroutine create_bc_mpi_datatype(nwstart, nwbc)
integer, dimension(:), allocatable recvrequest_c_sr
integer, dimension(-1:2^d &,-1:1^d &), target type_recv_srl_p1
integer, dimension(:), allocatable sendrequest_c_sr
integer, dimension(-1:2^d &,-1:1^d &), target type_send_srl_p2
integer, dimension(-1:2,-1:1) ixs_srl_
double precision, dimension(:), allocatable recvbuffer_srl
integer, dimension(:^d &,:^d &), pointer type_send_srl
integer, dimension(:,:), allocatable recvstatus_p
integer, dimension(^nd, 0:3) ixs_p_stg_
integer, dimension(^nd, 0:3^d &) sizes_p_recv_stg
integer, dimension(-1:1^d &, 0:3^d &), target type_recv_p_p2
integer, dimension(-1:1^d &, 0:3^d &), target type_recv_r_f
integer, dimension(:^d &,:^d &), pointer type_recv_r
integer, dimension(:^d &,:^d &), pointer type_send_r
integer, dimension(^nd,-1:1^d &) sizes_srl_recv_stg
double precision, dimension(:), allocatable recvbuffer_p
integer, dimension(^nd,-1:1) ixr_srl_stg_
integer, dimension(-1:1^d &, 0:3^d &), target type_recv_r_p1
double precision, dimension(:), allocatable sendbuffer_r
integer, dimension(-1:1^d &, 0:3^d &), target type_recv_p_f
integer, dimension(-1:1, 0:3) ixr_p_
integer, dimension(-1:1, 0:3) ixr_r_
integer, dimension(0:3^d &) sizes_p_recv_total
integer, dimension(:,:), allocatable recvstatus_c_sr
integer, dimension(0:3^d &) sizes_r_recv_total
integer, parameter npwbuf
integer, dimension(-1:1, 0:3) ixs_p_
integer, dimension(:), allocatable sendrequest_p
integer, dimension(^nd, 0:3^d &) sizes_r_recv_stg
subroutine getbc(time, qdt, psb, nwstart, nwbc, req_diag)
do update ghost cells of all blocks including physical boundaries
double precision, dimension(:), allocatable sendbuffer_srl
integer, dimension(-1:1^d &, 0:3^d &), target type_send_p_p2
integer, dimension(:,:), allocatable sendstatus_srl
integer, dimension(:,:), allocatable recvstatus_srl
integer, dimension(:), allocatable sendrequest_srl
integer, dimension(-1:2^d &,-1:1^d &), target type_send_srl_f
integer, dimension(:,:), allocatable recvstatus_r
integer, dimension(-1:1,-1:1) ixs_r_
integer, dimension(-1:1^d &,-1:1^d &), target type_send_r_p1
integer, dimension(0:3^d &) sizes_p_send_total
integer, dimension(:^d &,:^d &), pointer type_recv_p
integer, dimension(^nd, 0:3) ixr_p_stg_
integer, dimension(:^d &,:^d &), pointer type_send_p
integer, dimension(:,:), allocatable sendstatus_r
integer, dimension(:,:), allocatable sendstatus_p
integer, dimension(-1:1^d &) sizes_srl_recv_total
integer, dimension(-1:1^d &,-1:1^d &), target type_send_r_f
integer, dimension(:), allocatable recvrequest_p
integer, dimension(:), allocatable recvrequest_srl
integer, dimension(^nd,-1:1^d &) sizes_srl_send_stg
integer, dimension(:,:), allocatable recvstatus_c_p
This module contains definitions of global parameters and variables and some generic functions/subrou...
logical internalboundary
if there is an internal boundary
integer ixghi
Upper index of grid block arrays.
integer, dimension(3, 3) kr
Kronecker delta tensor.
integer, parameter ndim
Number of spatial dimensions for grid variables.
logical stagger_grid
True for using stagger grid.
logical, dimension(:), allocatable phyboundblock
True if a block has any physical boundary.
integer, dimension(:), allocatable, parameter d
logical ghost_copy
whether copy values instead of interpolation in ghost cells of finer blocks
integer ierrmpi
A global MPI error return code.
integer nghostcells
Number of ghost cells surrounding a grid.
This module defines the procedures of a physics module. It contains function pointers for the various...
Definition: mod_physics.t:4
procedure(sub_convert), pointer phys_to_primitive
Definition: mod_physics.t:59
logical phys_req_diagonal
Whether the physics routines require diagonal ghost cells, for example for computing a curl.
Definition: mod_physics.t:27
procedure(sub_convert), pointer phys_to_conserved
Definition: mod_physics.t:58
character(len=name_len) physics_type
String describing the physics type of the simulation.
Definition: mod_physics.t:19