MPI-AMRVAC  3.0
The MPI - Adaptive Mesh Refinement - Versatile Advection Code
mod_load_balance.t
Go to the documentation of this file.
2  implicit none
3 
4 contains
5  !> reallocate blocks into processors for load balance
6  subroutine load_balance
7  use mod_forest
10 
11  integer :: Morton_no, recv_igrid, recv_ipe, send_igrid, send_ipe, igrid, ipe
12  !> MPI recv send variables for AMR
13  integer :: itag, irecv, isend
14  integer, dimension(:), allocatable :: recvrequest, sendrequest
15  integer, dimension(:,:), allocatable :: recvstatus, sendstatus
16  !> MPI recv send variables for staggered-variable AMR
17  integer :: itag_stg
18  integer, dimension(:), allocatable :: recvrequest_stg, sendrequest_stg
19  integer, dimension(:,:), allocatable :: recvstatus_stg, sendstatus_stg
20 
21  integer, external :: getnode
22 
23  ! Jannis: for now, not using version for passive/active blocks
24  call get_morton_range()
25 
26  if (npe==1) then
28  return
29  end if
30 
31  irecv=0
32  isend=0
33  allocate(recvstatus(mpi_status_size,max_blocks),recvrequest(max_blocks), &
34  sendstatus(mpi_status_size,max_blocks),sendrequest(max_blocks))
35  recvrequest=mpi_request_null
36  sendrequest=mpi_request_null
37 
38  if(stagger_grid) then
39  allocate(recvstatus_stg(mpi_status_size,max_blocks*^nd),recvrequest_stg(max_blocks*^nd), &
40  sendstatus_stg(mpi_status_size,max_blocks*^nd),sendrequest_stg(max_blocks*^nd))
41  recvrequest_stg=mpi_request_null
42  sendrequest_stg=mpi_request_null
43  end if
44 
45  do ipe=0,npe-1; do morton_no=morton_start(ipe),morton_stop(ipe)
46  recv_ipe=ipe
47 
48  send_igrid=sfc(1,morton_no)
49  send_ipe=sfc(2,morton_no)
50 
51  if (recv_ipe/=send_ipe) then
52  ! get an igrid number for the new node in recv_ipe processor
53  recv_igrid=getnode(recv_ipe)
54  ! update node igrid and ipe on the tree
55  call change_ipe_tree_leaf(recv_igrid,recv_ipe,send_igrid,send_ipe)
56  ! receive physical data of the new node in recv_ipe processor
57  if (recv_ipe==mype) call lb_recv
58  ! send physical data of the old node in send_ipe processor
59  if (send_ipe==mype) call lb_send
60  end if
61  if (recv_ipe==mype) then
62  if (recv_ipe==send_ipe) then
63  sfc_to_igrid(morton_no)=send_igrid
64  else
65  sfc_to_igrid(morton_no)=recv_igrid
66  end if
67  end if
68  end do; end do
69 
70  if (irecv>0) then
71  call mpi_waitall(irecv,recvrequest,recvstatus,ierrmpi)
72  if(stagger_grid) call mpi_waitall(irecv,recvrequest_stg,recvstatus_stg,ierrmpi)
73  end if
74  if (isend>0) then
75  call mpi_waitall(isend,sendrequest,sendstatus,ierrmpi)
76  if(stagger_grid) call mpi_waitall(isend,sendrequest_stg,sendstatus_stg,ierrmpi)
77  end if
78 
79  deallocate(recvstatus,recvrequest,sendstatus,sendrequest)
80  if(stagger_grid) deallocate(recvstatus_stg,recvrequest_stg,sendstatus_stg,sendrequest_stg)
81 
82  ! post processing
83  do ipe=0,npe-1; do morton_no=morton_start(ipe),morton_stop(ipe)
84  recv_ipe=ipe
85 
86  send_igrid=sfc(1,morton_no)
87  send_ipe=sfc(2,morton_no)
88 
89  if (recv_ipe/=send_ipe) then
90  !if (send_ipe==mype) call dealloc_node(send_igrid)
91  call putnode(send_igrid,send_ipe)
92  end if
93  end do; end do
94  {#IFDEF EVOLVINGBOUNDARY
95  ! mark physical-boundary blocks on space-filling curve
96  do morton_no=morton_start(mype),morton_stop(mype)
97  igrid=sfc_to_igrid(morton_no)
98  if (phyboundblock(igrid)) sfc_phybound(morton_no)=1
99  end do
100  call mpi_allreduce(mpi_in_place,sfc_phybound,nleafs,mpi_integer,&
101  mpi_sum,icomm,ierrmpi)
102  }
103 
104  ! Update sfc array: igrid and ipe info in space filling curve
105  call amr_morton_order()
106 
107  contains
108 
109  subroutine lb_recv
110 
111  call alloc_node(recv_igrid)
112 
113  itag=recv_igrid
114  irecv=irecv+1
115  {#IFDEF EVOLVINGBOUNDARY
116  if (phyboundblock(recv_igrid)) then
117  call mpi_irecv(ps(recv_igrid)%w,1,type_block,send_ipe,itag, &
118  icomm,recvrequest(irecv),ierrmpi)
119  else
120  call mpi_irecv(ps(recv_igrid)%w,1,type_block_io,send_ipe,itag, &
121  icomm,recvrequest(irecv),ierrmpi)
122  end if
123  }{#IFNDEF EVOLVINGBOUNDARY
124  call mpi_irecv(ps(recv_igrid)%w,1,type_block_io,send_ipe,itag, &
125  icomm,recvrequest(irecv),ierrmpi)
126  }
127  if(stagger_grid) then
128  itag=recv_igrid+max_blocks
129  call mpi_irecv(ps(recv_igrid)%ws,1,type_block_io_stg,send_ipe,itag, &
130  icomm,recvrequest_stg(irecv),ierrmpi)
131  end if
132 
133  end subroutine lb_recv
134 
135  subroutine lb_send
136 
137  itag=recv_igrid
138  isend=isend+1
139  {#IFDEF EVOLVINGBOUNDARY
140  if (phyboundblock(send_igrid)) then
141  call mpi_isend(ps(send_igrid)%w,1,type_block,recv_ipe,itag, &
142  icomm,sendrequest(isend),ierrmpi)
143  else
144  call mpi_isend(ps(send_igrid)%w,1,type_block_io,recv_ipe,itag, &
145  icomm,sendrequest(isend),ierrmpi)
146  end if
147  }{#IFNDEF EVOLVINGBOUNDARY
148  call mpi_isend(ps(send_igrid)%w,1,type_block_io,recv_ipe,itag, &
149  icomm,sendrequest(isend),ierrmpi)
150  }
151  if(stagger_grid) then
152  itag=recv_igrid+max_blocks
153  call mpi_isend(ps(send_igrid)%ws,1,type_block_io_stg,recv_ipe,itag, &
154  icomm,sendrequest_stg(isend),ierrmpi)
155  end if
156 
157  end subroutine lb_send
158 
159  end subroutine load_balance
160 
161 end module mod_load_balance
subroutine alloc_node(igrid)
allocate arrays on igrid node
integer function getnode(ipe)
Get first available igrid on processor ipe.
subroutine putnode(igrid, ipe)
subroutine change_ipe_tree_leaf(recv_igrid, recv_ipe, send_igrid, send_ipe)
Definition: forest.t:218
subroutine lb_recv
subroutine lb_send
Module with basic grid data structures.
Definition: mod_forest.t:2
integer, dimension(:), allocatable, save sfc_to_igrid
Go from a Morton number to an igrid index (for a single processor)
Definition: mod_forest.t:53
integer, dimension(:), allocatable, save morton_start
First Morton number per processor.
Definition: mod_forest.t:62
integer, dimension(:), allocatable, save morton_stop
Last Morton number per processor.
Definition: mod_forest.t:65
integer, dimension(:,:), allocatable, save sfc
Array to go from a Morton number to an igrid and processor index. Sfc(1:3, MN) contains [igrid,...
Definition: mod_forest.t:43
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.
integer npe
The number of MPI tasks.
integer max_blocks
The maximum number of grid blocks in a processor.
subroutine load_balance
reallocate blocks into processors for load balance
subroutine get_morton_range
Set the Morton range for each processor.