41 integer :: iigrid, igrid
61 integer :: iigrid, igrid, i^
d, my_neighbor_type
62 integer :: iside, idim, ic^
d, inc^
d, ih^
d, icdim
64 integer :: idir,pi^
d, mi^
d, ph^
d, mh^
d, ipe_neighbor
65 integer :: nrecvs,nsends
67 integer :: nbuff_bc_recv_srl, nbuff_bc_send_srl, nbuff_bc_recv_r, nbuff_bc_send_r, nbuff_bc_recv_p, nbuff_bc_send_p
68 logical,
dimension(^ND) :: pole
75 nrecv_fc=0; nsend_fc=0
76 nbuff_bc_recv_srl=0; nbuff_bc_send_srl=0
77 nbuff_bc_recv_r=0; nbuff_bc_send_r=0
78 nbuff_bc_recv_p=0; nbuff_bc_send_p=0
81 do iigrid=1,igridstail; igrid=igrids(iigrid);
86 if (i^
d==0|.and.)
then
87 neighbor_type(0^
d&,igrid)=0
88 neighbor(1,0^
d&,igrid)=igrid
89 neighbor(2,0^
d&,igrid)=
mype
94 select case (my_neighbor_type)
96 case (neighbor_boundary)
97 neighbor(1,i^
d,igrid)=0
98 neighbor(2,i^
d,igrid)=-1
100 case (neighbor_coarse)
101 neighbor(1,i^
d,igrid)=my_neighbor%node%igrid
102 neighbor(2,i^
d,igrid)=my_neighbor%node%ipe
103 if (my_neighbor%node%ipe/=
mype)
then
104 ic^
d=1+modulo(tree%node%ig^
d-1,2);
105 if ({(i^
d==0.or.i^
d==2*ic^
d-3)|.and.})
then
115 case (neighbor_sibling)
116 neighbor(1,i^
d,igrid)=my_neighbor%node%igrid
117 neighbor(2,i^
d,igrid)=my_neighbor%node%ipe
118 if (my_neighbor%node%ipe/=
mype)
then
126 neighbor(1,i^
d,igrid)=0
127 neighbor(2,i^
d,igrid)=-1
130 {
do ic^db=1+int((1-i^db)/2),2-int((1+i^db)/2)
137 child%node => my_neighbor%node%child(ih^d)%node
138 neighbor_child(1,inc^d,igrid)=child%node%igrid
139 neighbor_child(2,inc^d,igrid)=child%node%ipe
140 if (child%node%ipe/=mype)
then
141 nrecv_bc_r=nrecv_bc_r+1
142 nsend_bc_p=nsend_bc_p+1
143 nbuff_bc_send_p=nbuff_bc_send_p+sizes_p_send_total(inc^d)
144 nbuff_bc_recv_r=nbuff_bc_recv_r+sizes_r_recv_total(inc^d)
150 if ({abs(i^d)+}==1)
then
155 select case (my_neighbor_type)
157 case (neighbor_coarse)
158 if (my_neighbor%node%ipe/=mype)
then
159 if (.not.pole(idim)) nsend_fc(idim)=nsend_fc(idim)+1
169 {
do ic^d=icdim,icdim^d%do ic^dd=1,2\}
170 child%node => my_neighbor%node%child(ic^dd)%node
171 if (child%node%ipe/=mype)
then
172 if (.not.pole(^d)) nrecv_fc(^d)=nrecv_fc(^d)+1
180 neighbor_pole(i^d,igrid)=0
181 if (my_neighbor_type>1)
then
184 neighbor_pole(i^d,igrid)=idim
190 neighbor_type(i^d,igrid)=my_neighbor_type
195 if(stagger_grid)
then
200 if ({abs(i^d)+}==1)
then
201 if (neighbor_pole(i^d,igrid)/=0) cycle
208 if (neighbor_type(i^d,igrid)==2)
then
210 pi^d=i^d+kr(idir,^d);
211 mi^d=i^d-kr(idir,^d);
212 ph^d=pi^d-kr(idim,^d)*(2*iside-3);
213 mh^d=mi^d-kr(idim,^d)*(2*iside-3);
215 if (neighbor_type(pi^d,igrid)==2.and.&
216 neighbor_type(ph^d,igrid)==2.and.&
217 mype/=neighbor(2,pi^d,igrid).and.&
218 neighbor_pole(pi^d,igrid)==0)
then
219 nsend_cc(idim) = nsend_cc(idim) + 1
222 if (neighbor_type(mi^d,igrid)==2.and.&
223 neighbor_type(mh^d,igrid)==2.and.&
224 mype/=neighbor(2,mi^d,igrid).and.&
225 neighbor_pole(mi^d,igrid)==0)
then
226 nsend_cc(idim) = nsend_cc(idim) + 1
231 if (neighbor_type(i^d,igrid)==3)
then
233 pi^d=i^d+kr(idir,^d);
234 mi^d=i^d-kr(idir,^d);
235 ph^d=pi^d-kr(idim,^d)*(2*iside-3);
236 mh^d=mi^d-kr(idim,^d)*(2*iside-3);
238 if (neighbor_type(pi^d,igrid)==4.and.&
239 neighbor_type(ph^d,igrid)==3.and.&
240 neighbor_pole(pi^d,igrid)==0)
then
242 {
do ic^db=1+int((1-pi^db)/2),2-int((1+pi^db)/2)
243 inc^db=2*pi^db+ic^db\}
244 if (mype.ne.neighbor_child(2,inc^d,igrid))
then
245 nrecv_cc(idim) = nrecv_cc(idim) + 1
250 if (neighbor_type(mi^d,igrid)==4.and.&
251 neighbor_type(mh^d,igrid)==3.and.&
252 neighbor_pole(mi^d,igrid)==0)
then
254 {
do ic^db=1+int((1-mi^db)/2),2-int((1+mi^db)/2)
255 inc^db=2*mi^db+ic^db\}
256 if (mype.ne.neighbor_child(2,inc^d,igrid))
then
257 nrecv_cc(idim) = nrecv_cc(idim) + 1
270 nrecvs=nrecv_bc_srl+nrecv_bc_r
271 if (
allocated(recvstatus_c_sr))
then
272 if(nrecvs >
size(recvrequest_c_sr))
then
273 deallocate(recvstatus_c_sr,recvrequest_c_sr)
274 allocate(recvstatus_c_sr(mpi_status_size,nrecvs),recvrequest_c_sr(nrecvs))
277 allocate(recvstatus_c_sr(mpi_status_size,nrecvs),recvrequest_c_sr(nrecvs))
279 recvrequest_c_sr=mpi_request_null
282 nsends=nsend_bc_srl+nsend_bc_r
283 if (
allocated(sendstatus_c_sr))
then
284 if(nsends >
size(sendrequest_c_sr))
then
285 deallocate(sendstatus_c_sr,sendrequest_c_sr)
286 allocate(sendstatus_c_sr(mpi_status_size,nsends),sendrequest_c_sr(nsends))
289 allocate(sendstatus_c_sr(mpi_status_size,nsends),sendrequest_c_sr(nsends))
291 sendrequest_c_sr=mpi_request_null
294 if (
allocated(recvstatus_c_p))
then
295 if(nrecv_bc_p >
size(recvrequest_c_p))
then
296 deallocate(recvstatus_c_p,recvrequest_c_p)
297 allocate(recvstatus_c_p(mpi_status_size,nrecv_bc_p),recvrequest_c_p(nrecv_bc_p))
300 allocate(recvstatus_c_p(mpi_status_size,nrecv_bc_p),recvrequest_c_p(nrecv_bc_p))
302 recvrequest_c_p=mpi_request_null
305 if (
allocated(sendstatus_c_p))
then
306 if(nsend_bc_p >
size(sendrequest_c_p))
then
307 deallocate(sendstatus_c_p,sendrequest_c_p)
308 allocate(sendstatus_c_p(mpi_status_size,nsend_bc_p),sendrequest_c_p(nsend_bc_p))
311 allocate(sendstatus_c_p(mpi_status_size,nsend_bc_p),sendrequest_c_p(nsend_bc_p))
313 sendrequest_c_p=mpi_request_null
315 if(stagger_grid)
then
317 if (
allocated(recvbuffer_srl))
then
318 if (nbuff_bc_recv_srl >
size(recvbuffer_srl))
then
319 deallocate(recvbuffer_srl)
320 allocate(recvbuffer_srl(nbuff_bc_recv_srl))
323 allocate(recvbuffer_srl(nbuff_bc_recv_srl))
325 if (
allocated(recvstatus_srl))
then
326 if(nrecv_bc_srl>
size(recvrequest_srl))
then
327 deallocate(recvstatus_srl,recvrequest_srl)
328 allocate(recvstatus_srl(mpi_status_size,nrecv_bc_srl),recvrequest_srl(nrecv_bc_srl))
331 allocate(recvstatus_srl(mpi_status_size,nrecv_bc_srl),recvrequest_srl(nrecv_bc_srl))
333 recvrequest_srl=mpi_request_null
336 if (
allocated(sendbuffer_srl))
then
337 if (nbuff_bc_send_srl >
size(sendbuffer_srl))
then
338 deallocate(sendbuffer_srl)
339 allocate(sendbuffer_srl(nbuff_bc_send_srl))
342 allocate(sendbuffer_srl(nbuff_bc_send_srl))
344 if (
allocated(sendstatus_srl))
then
345 if(nsend_bc_srl>
size(sendrequest_srl))
then
346 deallocate(sendstatus_srl,sendrequest_srl)
347 allocate(sendstatus_srl(mpi_status_size,nsend_bc_srl),sendrequest_srl(nsend_bc_srl))
350 allocate(sendstatus_srl(mpi_status_size,nsend_bc_srl),sendrequest_srl(nsend_bc_srl))
352 sendrequest_srl=mpi_request_null
355 if (
allocated(recvbuffer_r))
then
356 if (nbuff_bc_recv_r >
size(recvbuffer_r))
then
357 deallocate(recvbuffer_r)
358 allocate(recvbuffer_r(nbuff_bc_recv_r))
361 allocate(recvbuffer_r(nbuff_bc_recv_r))
363 if (
allocated(recvstatus_r))
then
364 if(nrecv_bc_r>
size(recvrequest_r))
then
365 deallocate(recvstatus_r,recvrequest_r)
366 allocate(recvstatus_r(mpi_status_size,nrecv_bc_r),recvrequest_r(nrecv_bc_r))
369 allocate(recvstatus_r(mpi_status_size,nrecv_bc_r),recvrequest_r(nrecv_bc_r))
371 recvrequest_r=mpi_request_null
374 if (
allocated(sendbuffer_r))
then
375 if (nbuff_bc_send_r >
size(sendbuffer_r))
then
376 deallocate(sendbuffer_r)
377 allocate(sendbuffer_r(nbuff_bc_send_r))
380 allocate(sendbuffer_r(nbuff_bc_send_r))
382 if (
allocated(sendstatus_r))
then
383 if(nsend_bc_r>
size(sendrequest_r))
then
384 deallocate(sendstatus_r,sendrequest_r)
385 allocate(sendstatus_r(mpi_status_size,nsend_bc_r),sendrequest_r(nsend_bc_r))
388 allocate(sendstatus_r(mpi_status_size,nsend_bc_r),sendrequest_r(nsend_bc_r))
390 sendrequest_r=mpi_request_null
393 if (
allocated(recvbuffer_p))
then
394 if (nbuff_bc_recv_p >
size(recvbuffer_p))
then
395 deallocate(recvbuffer_p)
396 allocate(recvbuffer_p(nbuff_bc_recv_p))
399 allocate(recvbuffer_p(nbuff_bc_recv_p))
401 if (
allocated(recvstatus_p))
then
402 if(nrecv_bc_p>
size(recvrequest_p))
then
403 deallocate(recvstatus_p,recvrequest_p)
404 allocate(recvstatus_p(mpi_status_size,nrecv_bc_p),recvrequest_p(nrecv_bc_p))
407 allocate(recvstatus_p(mpi_status_size,nrecv_bc_p),recvrequest_p(nrecv_bc_p))
409 recvrequest_p=mpi_request_null
412 if (
allocated(sendbuffer_p))
then
413 if (nbuff_bc_send_p >
size(sendbuffer_p))
then
414 deallocate(sendbuffer_p)
415 allocate(sendbuffer_p(nbuff_bc_send_p))
418 allocate(sendbuffer_p(nbuff_bc_send_p))
420 if (
allocated(sendstatus_p))
then
421 if(nsend_bc_p>
size(sendrequest_p))
then
422 deallocate(sendstatus_p,sendrequest_p)
423 allocate(sendstatus_p(mpi_status_size,nsend_bc_p),sendrequest_p(nsend_bc_p))
426 allocate(sendstatus_p(mpi_status_size,nsend_bc_p),sendrequest_p(nsend_bc_p))
428 sendrequest_p=mpi_request_null
subroutine, public find_neighbor(my_neighbor, my_neighbor_type, tree, iD, pole)
find neighors of all blocks
Module with basic grid data structures.
logical, dimension(:,:), allocatable, save igrid_inuse
type(tree_node_ptr), dimension(:), allocatable, save level_tail
The tail pointer of the linked list per refinement level.
type(tree_node_ptr), dimension(:,:), allocatable, save igrid_to_node
Array to go from an [igrid, ipe] index to a node pointer.
subroutine, public get_level_range
subroutine, public build_connectivity
subroutine, public getigrids
update ghost cells of all blocks including physical boundaries
integer, dimension(-1:1^d &) sizes_r_send_total
integer, dimension(-1:1^d &) sizes_srl_send_total
integer, dimension(0:3^d &) sizes_p_recv_total
integer, dimension(-1:1^d &) sizes_srl_recv_total
This module contains definitions of global parameters and variables and some generic functions/subrou...
logical stagger_grid
True for using stagger grid.
integer mype
The rank of the current MPI task.
double precision, dimension(:), allocatable, parameter d
integer refine_max_level
Maximal number of AMR levels.
integer max_blocks
The maximum number of grid blocks in a processor.