MPI-AMRVAC  3.1
The MPI - Adaptive Mesh Refinement - Versatile Advection Code (development version)
mod_functions_connectivity.t
Go to the documentation of this file.
2 
3 
4  implicit none
5  private
6 
7  public :: get_level_range
8  public :: getigrids
9  public :: build_connectivity
10  contains
11 
12  subroutine get_level_range
13  use mod_forest
15 
16  integer :: level
17 
18  ! determine new finest level
19  do level=refine_max_level,1,-1
20  if (associated(level_tail(level)%node)) then
21  levmax=level
22  exit
23  end if
24  end do
25 
26  ! determine coarsest level
27  do level=1,levmax
28  if (associated(level_tail(level)%node)) then
29  levmin=level
30  exit
31  end if
32  end do
33 
34  end subroutine get_level_range
35 
36  subroutine getigrids
37  use mod_forest
39  implicit none
40 
41  integer :: iigrid, igrid
42 
43  iigrid=0
44  do igrid=1,max_blocks
45  if (igrid_inuse(igrid,mype)) then
46  iigrid=iigrid+1
47  igrids(iigrid)=igrid
48  end if
49  end do
50 
51  igridstail=iigrid
52 
53  end subroutine getigrids
54 
55  subroutine build_connectivity
56  use mod_forest
60 
61  integer :: iigrid, igrid, i^d, my_neighbor_type
62  integer :: iside, idim, ic^d, inc^d, ih^d, icdim
63  ! Variables to detect special corners for stagger grid
64  integer :: idir,pi^d, mi^d, ph^d, mh^d, ipe_neighbor
65  integer :: nrecvs,nsends
66  ! total size of buffer arrays
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
69  logical :: nopole
70  type(tree_node_ptr) :: tree, my_neighbor, child
71 
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
79  if(stagger_grid) nrecv_cc=0; nsend_cc=0
80 
81  do iigrid=1,igridstail; igrid=igrids(iigrid);
82  tree%node => igrid_to_node(igrid,mype)%node
83 
84  {do i^db=-1,1\}
85  ! skip the grid itself
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
90  else
91  call find_neighbor(my_neighbor,my_neighbor_type,tree,i^d,pole)
92  nopole=.not.any(pole)
93 
94  select case (my_neighbor_type)
95  ! adjacent to physical boundary
96  case (neighbor_boundary)
97  neighbor(1,i^d,igrid)=0
98  neighbor(2,i^d,igrid)=-1
99  ! fine-coarse transition
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
108  nbuff_bc_send_r=nbuff_bc_send_r+sizes_r_send_total(i^d)
109  ! This is the local index of the prolonged ghost zone
110  inc^d=ic^d+i^d;
111  nbuff_bc_recv_p=nbuff_bc_recv_p+sizes_p_recv_total(inc^d)
112  end if
113  end if
114  ! same refinement level
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
121  nbuff_bc_send_srl=nbuff_bc_send_srl+sizes_srl_send_total(i^d)
122  nbuff_bc_recv_srl=nbuff_bc_recv_srl+sizes_srl_recv_total(i^d)
123  end if
124  ! coarse-fine transition
125  case (neighbor_fine)
126  neighbor(1,i^d,igrid)=0
127  neighbor(2,i^d,igrid)=-1
128  ! Loop over the local indices of children ic^D
129  ! and calculate local indices of ghost zone inc^D.
130  {do ic^db=1+int((1-i^db)/2),2-int((1+i^db)/2)
131  inc^db=2*i^db+ic^db
132  if (pole(^db)) then
133  ih^db=3-ic^db
134  else
135  ih^db=ic^db
136  end if\}
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)
145  end if
146  {end do\}
147  end select
148 
149  ! flux fix for conservation only for pure directional shifts
150  if ({abs(i^d)+}==1) then
151  {if (i^d/=0) then
152  idim=^d
153  iside=int((i^d+3)/2)
154  end if\}
155  select case (my_neighbor_type)
156  ! only across fine-coarse or coarse-fine boundaries
157  case (neighbor_coarse)
158  if (my_neighbor%node%ipe/=mype) then
159  if (.not.pole(idim)) nsend_fc(idim)=nsend_fc(idim)+1
160  end if
161  case (neighbor_fine)
162  if (pole(idim)) then
163  icdim=iside
164  else
165  icdim=3-iside
166  end if
167  select case (idim)
168  {case (^d)
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
173  end if
174  {end do\} \}
175  end select
176  end select
177  end if
178 
179  if (phi_ > 0) then
180  neighbor_pole(i^d,igrid)=0
181  if (my_neighbor_type>1) then
182  do idim=1,^nd
183  if (pole(idim)) then
184  neighbor_pole(i^d,igrid)=idim
185  exit ! there can only be one pole between two meshes
186  end if
187  end do
188  end if
189  end if
190  neighbor_type(i^d,igrid)=my_neighbor_type
191 
192  end if
193  {end do\}
194 
195  if(stagger_grid) then
196  !Now all the neighbour information is known.
197  !Check if there are special corners that need to be communicated
198  !To determine whether to send/receive, we must check three neighbours
199  {do i^db=-1,1\}
200  if ({abs(i^d)+}==1) then
201  if (neighbor_pole(i^d,igrid)/=0) cycle
202  ! Assign value to idim and iside
203  {if (i^d/=0) then
204  idim=^d
205  iside=int((i^d+3)/2)
206  end if\}
207  ! Fine block surrounded by coarse blocks
208  if (neighbor_type(i^d,igrid)==2) then
209  do idir=idim+1,ndim
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);
214 
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
220  end if
221 
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
227  end if
228  end do
229  end if
230  ! Coarse block diagonal to fine block(s)
231  if (neighbor_type(i^d,igrid)==3) then
232  do idir=idim+1,ndim
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);
237 
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
241  ! Loop on children (several in 3D)
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
246  end if
247  {end do\}
248  end if
249 
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
253  ! Loop on children (several in 3D)
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
258  end if
259  {end do\}
260  end if
261  end do
262  end if
263  end if
264  {end do\}
265  end if
266 
267  end do
268 
269  ! allocate space for mpi recieve for siblings and restrict ghost cell filling
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))
275  end if
276  else
277  allocate(recvstatus_c_sr(mpi_status_size,nrecvs),recvrequest_c_sr(nrecvs))
278  end if
279  recvrequest_c_sr=mpi_request_null
280 
281  ! allocate space for mpi send for siblings and restrict ghost cell filling
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))
287  end if
288  else
289  allocate(sendstatus_c_sr(mpi_status_size,nsends),sendrequest_c_sr(nsends))
290  end if
291  sendrequest_c_sr=mpi_request_null
292 
293  ! allocate space for mpi receive for prolongation ghost cell filling
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))
298  end if
299  else
300  allocate(recvstatus_c_p(mpi_status_size,nrecv_bc_p),recvrequest_c_p(nrecv_bc_p))
301  end if
302  recvrequest_c_p=mpi_request_null
303 
304  ! allocate space for mpi send for prolongation ghost cell filling
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))
309  end if
310  else
311  allocate(sendstatus_c_p(mpi_status_size,nsend_bc_p),sendrequest_c_p(nsend_bc_p))
312  end if
313  sendrequest_c_p=mpi_request_null
314 
315  if(stagger_grid) then
316  ! allocate space for receive buffer for siblings ghost cell filling
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))
321  end if
322  else
323  allocate(recvbuffer_srl(nbuff_bc_recv_srl))
324  end if
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))
329  end if
330  else
331  allocate(recvstatus_srl(mpi_status_size,nrecv_bc_srl),recvrequest_srl(nrecv_bc_srl))
332  end if
333  recvrequest_srl=mpi_request_null
334 
335  ! allocate space for send buffer for siblings ghost cell filling
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))
340  end if
341  else
342  allocate(sendbuffer_srl(nbuff_bc_send_srl))
343  end if
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))
348  end if
349  else
350  allocate(sendstatus_srl(mpi_status_size,nsend_bc_srl),sendrequest_srl(nsend_bc_srl))
351  end if
352  sendrequest_srl=mpi_request_null
353 
354  ! allocate space for recieve buffer for restrict ghost cell filling
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))
359  end if
360  else
361  allocate(recvbuffer_r(nbuff_bc_recv_r))
362  end if
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))
367  end if
368  else
369  allocate(recvstatus_r(mpi_status_size,nrecv_bc_r),recvrequest_r(nrecv_bc_r))
370  end if
371  recvrequest_r=mpi_request_null
372 
373  ! allocate space for send buffer for restrict ghost cell filling
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))
378  end if
379  else
380  allocate(sendbuffer_r(nbuff_bc_send_r))
381  end if
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))
386  end if
387  else
388  allocate(sendstatus_r(mpi_status_size,nsend_bc_r),sendrequest_r(nsend_bc_r))
389  end if
390  sendrequest_r=mpi_request_null
391 
392  ! allocate space for recieve buffer for prolong ghost cell filling
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))
397  end if
398  else
399  allocate(recvbuffer_p(nbuff_bc_recv_p))
400  end if
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))
405  end if
406  else
407  allocate(recvstatus_p(mpi_status_size,nrecv_bc_p),recvrequest_p(nrecv_bc_p))
408  end if
409  recvrequest_p=mpi_request_null
410 
411  ! allocate space for send buffer for restrict ghost cell filling
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))
416  end if
417  else
418  allocate(sendbuffer_p(nbuff_bc_send_p))
419  end if
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))
424  end if
425  else
426  allocate(sendstatus_p(mpi_status_size,nsend_bc_p),sendrequest_p(nsend_bc_p))
427  end if
428  sendrequest_p=mpi_request_null
429  end if
430 
431  end subroutine build_connectivity
432 
subroutine, public find_neighbor(my_neighbor, my_neighbor_type, tree, iD, pole)
find neighors of all blocks
Module with basic grid data structures.
Definition: mod_forest.t:2
logical, dimension(:,:), allocatable, save igrid_inuse
Definition: mod_forest.t:70
type(tree_node_ptr), dimension(:), allocatable, save level_tail
The tail pointer of the linked list per refinement level.
Definition: mod_forest.t:38
type(tree_node_ptr), dimension(:,:), allocatable, save igrid_to_node
Array to go from an [igrid, ipe] index to a node pointer.
Definition: mod_forest.t:32
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.
Pointer to a tree_node.
Definition: mod_forest.t:6