1 | ! This module contains the requests for transfert_mpi |
---|
2 | module transfert_request_mod |
---|
3 | use abort_mod |
---|
4 | use field_mod, only : t_field, field_T, field_U |
---|
5 | use domain_mod, only : ndomain, ndomain_glo, domain, domain_glo, domloc_glo_ind, domglo_rank, domglo_loc_ind |
---|
6 | implicit none |
---|
7 | private |
---|
8 | |
---|
9 | ! Points for t_request : list of points coordinates in a domain |
---|
10 | type t_points |
---|
11 | integer, allocatable :: i(:), j(:), edge(:) |
---|
12 | integer :: npoints |
---|
13 | end type |
---|
14 | |
---|
15 | ! Describes how to exchange data from/to one domain to/from other domains |
---|
16 | type t_request |
---|
17 | integer :: field_type |
---|
18 | type(t_points), allocatable :: points_BtoH(:) ! req(local_ind)%points_BtoH(remote_ind) are the points of |
---|
19 | ! domain local_ind to unpack from mpi buffer recieved from remote_ind |
---|
20 | type(t_points), allocatable :: points_HtoB(:) ! req(local_ind)%points_HtoB(remote_ind) are the points of |
---|
21 | ! domain local_ind to pack to mpi buffer to send to remote_ind |
---|
22 | logical :: vector = .false. ! Is a vector request (true if sign change needed) |
---|
23 | end type t_request |
---|
24 | |
---|
25 | type(t_request),save,pointer :: req_i1(:) ! Halos for T fields |
---|
26 | type(t_request),save,pointer :: req_e1_scal(:) ! Halos for U fields |
---|
27 | type(t_request),save,pointer :: req_e1_vect(:) ! Halos for U fields (with sign change) |
---|
28 | type(t_request),save,pointer :: req_z1_scal(:) => null() ! DO NOT USE : not initialized |
---|
29 | |
---|
30 | type(t_request),save,pointer :: req_i0(:) ! Duplicated cells for T fields |
---|
31 | type(t_request),save,pointer :: req_e0_scal(:) |
---|
32 | type(t_request),save,pointer :: req_e0_vect(:) |
---|
33 | |
---|
34 | public :: t_points, t_request, & |
---|
35 | req_i1, req_e1_scal, req_e1_vect, & |
---|
36 | req_i0, req_e0_scal, req_e0_vect, & |
---|
37 | req_z1_scal, & |
---|
38 | init_all_requests |
---|
39 | contains |
---|
40 | |
---|
41 | ! Initialize requests |
---|
42 | subroutine init_all_requests |
---|
43 | use dimensions, only : swap_dimensions, ii_begin, ii_end, jj_begin, jj_end |
---|
44 | use metric, only : right, rdown, ldown, left, lup, rup |
---|
45 | integer :: ind, i, j |
---|
46 | |
---|
47 | ! Create requests req_* see corresponding module variables for description |
---|
48 | call request_init(req_i0, field_T) |
---|
49 | DO ind=1,ndomain |
---|
50 | CALL swap_dimensions(ind) |
---|
51 | DO i=ii_begin,ii_end |
---|
52 | CALL request_add_point(req_i0,ind,i,jj_begin) |
---|
53 | ENDDO |
---|
54 | DO j=jj_begin,jj_end |
---|
55 | CALL request_add_point(req_i0,ind,ii_begin,j) |
---|
56 | CALL request_add_point(req_i0,ind,ii_end,j) |
---|
57 | ENDDO |
---|
58 | DO i=ii_begin,ii_end |
---|
59 | CALL request_add_point(req_i0,ind,i,jj_end) |
---|
60 | ENDDO |
---|
61 | ENDDO |
---|
62 | call request_exchange(req_i0) |
---|
63 | |
---|
64 | call request_init(req_i1, field_T) |
---|
65 | do ind=1,ndomain |
---|
66 | call swap_dimensions(ind) |
---|
67 | do i=ii_begin,ii_end+1 |
---|
68 | call request_add_point(req_i1,ind,i,jj_begin-1) |
---|
69 | enddo |
---|
70 | do j=jj_begin,jj_end |
---|
71 | call request_add_point(req_i1,ind,ii_begin-1,j) |
---|
72 | call request_add_point(req_i1,ind,ii_end+1,j) |
---|
73 | enddo |
---|
74 | call request_add_point(req_i1,ind,ii_begin-1,jj_end+1) |
---|
75 | do i=ii_begin,ii_end |
---|
76 | call request_add_point(req_i1,ind,i,jj_end+1) |
---|
77 | enddo |
---|
78 | enddo |
---|
79 | call request_exchange(req_i1) |
---|
80 | |
---|
81 | call request_init(req_e0_scal, field_U) |
---|
82 | DO ind=1,ndomain |
---|
83 | CALL swap_dimensions(ind) |
---|
84 | DO i=ii_begin+1,ii_end-1 |
---|
85 | CALL request_add_point(req_e0_scal,ind,i,jj_begin,right) |
---|
86 | CALL request_add_point(req_e0_scal,ind,i,jj_end,right) |
---|
87 | ENDDO |
---|
88 | DO j=jj_begin+1,jj_end-1 |
---|
89 | CALL request_add_point(req_e0_scal,ind,ii_begin,j,rup) |
---|
90 | CALL request_add_point(req_e0_scal,ind,ii_end,j,rup) |
---|
91 | ENDDO |
---|
92 | CALL request_add_point(req_e0_scal,ind,ii_begin+1,jj_begin,left) |
---|
93 | CALL request_add_point(req_e0_scal,ind,ii_begin,jj_begin+1,ldown) |
---|
94 | CALL request_add_point(req_e0_scal,ind,ii_begin+1,jj_end,left) |
---|
95 | CALL request_add_point(req_e0_scal,ind,ii_end,jj_begin+1,ldown) |
---|
96 | ENDDO |
---|
97 | call request_exchange(req_e0_scal) |
---|
98 | |
---|
99 | ! req_e0_vect is the same request than req_e0_scal but with the vector flag to true |
---|
100 | allocate( req_e0_vect(ndomain) ) |
---|
101 | req_e0_vect = req_e0_scal |
---|
102 | req_e0_vect%vector = .true. |
---|
103 | |
---|
104 | call request_init(req_e1_scal, field_U) |
---|
105 | DO ind=1,ndomain |
---|
106 | CALL swap_dimensions(ind) |
---|
107 | DO i=ii_begin,ii_end |
---|
108 | CALL request_add_point(req_e1_scal,ind,i,jj_begin-1,rup) |
---|
109 | CALL request_add_point(req_e1_scal,ind,i+1,jj_begin-1,lup) |
---|
110 | ENDDO |
---|
111 | DO j=jj_begin,jj_end |
---|
112 | CALL request_add_point(req_e1_scal,ind,ii_end+1,j,left) |
---|
113 | ENDDO |
---|
114 | DO j=jj_begin,jj_end |
---|
115 | CALL request_add_point(req_e1_scal,ind,ii_end+1,j-1,lup) |
---|
116 | ENDDO |
---|
117 | DO i=ii_begin,ii_end |
---|
118 | CALL request_add_point(req_e1_scal,ind,i,jj_end+1,ldown) |
---|
119 | CALL request_add_point(req_e1_scal,ind,i-1,jj_end+1,rdown) |
---|
120 | ENDDO |
---|
121 | DO j=jj_begin,jj_end |
---|
122 | CALL request_add_point(req_e1_scal,ind,ii_begin-1,j,right) |
---|
123 | ENDDO |
---|
124 | DO j=jj_begin,jj_end |
---|
125 | CALL request_add_point(req_e1_scal,ind,ii_begin-1,j+1,rdown) |
---|
126 | ENDDO |
---|
127 | ENDDO |
---|
128 | call request_exchange(req_e1_scal) |
---|
129 | |
---|
130 | allocate( req_e1_vect(ndomain) ) |
---|
131 | req_e1_vect = req_e1_scal |
---|
132 | req_e1_vect%vector = .true. |
---|
133 | contains |
---|
134 | |
---|
135 | ! ---- Points ---- |
---|
136 | ! Initialize (allocate) points |
---|
137 | subroutine points_init(points, has_edges) |
---|
138 | type(t_points):: points |
---|
139 | logical :: has_edges |
---|
140 | integer, parameter :: INITIAL_ALLOC_SIZE = 10 |
---|
141 | |
---|
142 | allocate(points%i(INITIAL_ALLOC_SIZE)) |
---|
143 | allocate(points%j(INITIAL_ALLOC_SIZE)) |
---|
144 | if(has_edges) allocate(points%edge(INITIAL_ALLOC_SIZE)) |
---|
145 | points%npoints = 0 |
---|
146 | end subroutine |
---|
147 | |
---|
148 | ! Append point (i, j [,edge]) to point list |
---|
149 | subroutine points_append(points, i, j, edge) |
---|
150 | type(t_points), intent(inout):: points |
---|
151 | integer, intent(in) :: i, j |
---|
152 | integer, intent(in), optional :: edge |
---|
153 | integer :: begin_size |
---|
154 | |
---|
155 | begin_size=points%npoints |
---|
156 | call array_append( points%i, begin_size, i ) |
---|
157 | begin_size=points%npoints |
---|
158 | call array_append( points%j, begin_size, j ) |
---|
159 | if(present(edge)) then |
---|
160 | begin_size=points%npoints |
---|
161 | call array_append( points%edge, begin_size, edge ) |
---|
162 | end if |
---|
163 | |
---|
164 | points%npoints = points%npoints + 1 |
---|
165 | end subroutine |
---|
166 | |
---|
167 | ! Add element to array, and reallocate if necessary |
---|
168 | subroutine array_append( a, a_size, elt ) |
---|
169 | integer, allocatable, intent(inout) :: a(:) |
---|
170 | integer, intent(inout) :: a_size |
---|
171 | integer, intent(in) :: elt |
---|
172 | integer, allocatable :: a_tmp(:) |
---|
173 | integer, parameter :: GROW_FACTOR = 2 |
---|
174 | |
---|
175 | if( size( a ) <= a_size ) then |
---|
176 | allocate( a_tmp ( a_size * GROW_FACTOR ) ) |
---|
177 | a_tmp(1:a_size) = a(1:a_size) |
---|
178 | call move_alloc(a_tmp, a) |
---|
179 | end if |
---|
180 | a_size = a_size + 1 |
---|
181 | a(a_size) = elt; |
---|
182 | end subroutine |
---|
183 | |
---|
184 | ! Sort and remove duplicates in points |
---|
185 | ! dimensions%iim must be set to the right value with swap_dimensions before |
---|
186 | ! After this call, size(points%i) == points%npoints |
---|
187 | subroutine points_sort( points ) |
---|
188 | use dimensions, only : iim, jjm |
---|
189 | type(t_points):: points |
---|
190 | integer :: displs(points%npoints) |
---|
191 | integer :: k, i, last, i_min, tmp |
---|
192 | |
---|
193 | displs = (points%i-1) + (points%j-1)*iim |
---|
194 | if(allocated(points%edge)) displs = displs + (points%edge-1)*iim*jjm |
---|
195 | |
---|
196 | last=0 |
---|
197 | do i=1, points%npoints |
---|
198 | i_min = minloc(displs(i:),1)+i-1 |
---|
199 | tmp = displs(i_min) |
---|
200 | displs(i_min) = displs(i) |
---|
201 | if(last==0 .or. displs(i_min) /= displs(last)) last=last+1 |
---|
202 | displs(last) = tmp |
---|
203 | end do |
---|
204 | |
---|
205 | points%npoints = last |
---|
206 | deallocate(points%i); allocate(points%i(last)) |
---|
207 | deallocate(points%j); allocate(points%j(last)) |
---|
208 | do k=1, last |
---|
209 | points%i(k) = mod(displs(k), iim)+1 |
---|
210 | points%j(k) = mod((displs(k)-(points%i(k)-1))/iim, jjm)+1 |
---|
211 | end do |
---|
212 | if(allocated(points%edge)) then |
---|
213 | deallocate(points%edge); allocate(points%edge(last)) |
---|
214 | points%edge(:) = displs(1:last)/(iim*jjm) + 1 |
---|
215 | end if |
---|
216 | |
---|
217 | if( size(points%i) /= points%npoints ) call dynamico_abort( "Storage too big for points" ) |
---|
218 | end subroutine |
---|
219 | |
---|
220 | ! ---- Requests ---- |
---|
221 | ! Initialize request |
---|
222 | subroutine request_init( request, field_type ) |
---|
223 | type(t_request), pointer, intent(out) :: request(:) |
---|
224 | integer, intent(in) :: field_type |
---|
225 | integer :: ind, remote_ind_glo |
---|
226 | |
---|
227 | allocate( request(ndomain) ) |
---|
228 | do ind=1, ndomain |
---|
229 | request(ind)%field_type = field_type |
---|
230 | allocate(request(ind)%points_BtoH(ndomain_glo)) |
---|
231 | allocate(request(ind)%points_HtoB(ndomain_glo)) |
---|
232 | do remote_ind_glo = 1, ndomain_glo |
---|
233 | call points_init(request(ind)%points_BtoH(remote_ind_glo), field_type == field_U) |
---|
234 | !call points_init(request(ind)%points_HtoB(remote_ind_glo)) ! Initialized before MPI communication in request_exchange |
---|
235 | end do |
---|
236 | end do |
---|
237 | end subroutine |
---|
238 | |
---|
239 | ! Append local point (ind_glo,i,j) to the point list of domain ind, where ind_glo is the domain owning the point |
---|
240 | subroutine request_add_point(req, ind, i, j, edge) |
---|
241 | type(t_request), intent(inout) :: req(:) |
---|
242 | integer :: ind, i, j |
---|
243 | integer, optional :: edge |
---|
244 | |
---|
245 | if( present(edge) ) then ! U field |
---|
246 | if(domain(ind)%edge_assign_domain(edge-1,i,j) /= domloc_glo_ind(ind) ) & ! Add only if remote domain != local domain |
---|
247 | call points_append(req(ind)%points_BtoH(domain(ind)%edge_assign_domain(edge-1,i,j)), i, j, edge ) |
---|
248 | else ! T field |
---|
249 | if(domain(ind)%assign_domain(i,j) /= domloc_glo_ind(ind) ) & |
---|
250 | call points_append(req(ind)%points_BtoH(domain(ind)%assign_domain(i,j)), i, j ) |
---|
251 | endif |
---|
252 | end subroutine |
---|
253 | |
---|
254 | ! Each domain send the cell index it needs from other domains |
---|
255 | subroutine request_exchange(req) |
---|
256 | use mpi_mod |
---|
257 | use dimensions, only : swap_dimensions |
---|
258 | use mpipara, only : comm_icosa |
---|
259 | type(t_request), target, intent(in) :: req(:) |
---|
260 | type(t_points) :: points_send(ndomain,ndomain_glo) |
---|
261 | type(t_points), pointer :: points |
---|
262 | |
---|
263 | integer :: ind_loc, remote_ind_glo, k |
---|
264 | integer :: i, j, edge |
---|
265 | integer :: reqs(ndomain, ndomain_glo, 2), reqs2(ndomain, ndomain_glo, 2), reqs3(ndomain, ndomain_glo, 3), ierr |
---|
266 | |
---|
267 | ! Transform position in local domain to position in remote domain pos : (i,j) -> pos_send : (assign_i(i,j), assign_j(i,j)) |
---|
268 | do ind_loc = 1, ndomain |
---|
269 | call swap_dimensions(ind_loc) |
---|
270 | do remote_ind_glo = 1, ndomain_glo |
---|
271 | points => req(ind_loc)%points_BtoH(remote_ind_glo) |
---|
272 | call points_sort(points) |
---|
273 | call points_init( points_send(ind_loc,remote_ind_glo), req(ind_loc)%field_type == field_U ) |
---|
274 | do k = 1, points%npoints |
---|
275 | i = points%i(k) |
---|
276 | j = points%j(k) |
---|
277 | if(req(ind_loc)%field_type == field_T) then |
---|
278 | call points_append(points_send(ind_loc,remote_ind_glo), domain(ind_loc)%assign_i(i,j), & |
---|
279 | domain(ind_loc)%assign_j(i,j)) |
---|
280 | else |
---|
281 | edge = points%edge(k) - 1 ! WARNING : domain%edge_assign_* use index in 0:5 instead of 1:6 everywhere else |
---|
282 | call points_append(points_send(ind_loc,remote_ind_glo), domain(ind_loc)%edge_assign_i(edge,i,j), & |
---|
283 | domain(ind_loc)%edge_assign_j(edge,i,j), & |
---|
284 | domain(ind_loc)%edge_assign_pos(edge,i,j)+1) |
---|
285 | end if |
---|
286 | end do |
---|
287 | end do |
---|
288 | end do |
---|
289 | |
---|
290 | ! Exchange sizes |
---|
291 | do ind_loc = 1, ndomain |
---|
292 | do remote_ind_glo = 1, ndomain_glo |
---|
293 | call MPI_Isend( points_send(ind_loc,remote_ind_glo)%npoints, 1, MPI_INT, domglo_rank(remote_ind_glo), & |
---|
294 | remote_ind_glo+ndomain_glo*domloc_glo_ind(ind_loc), comm_icosa, reqs(ind_loc,remote_ind_glo,1), ierr ) |
---|
295 | call MPI_Irecv( req(ind_loc)%points_HtoB(remote_ind_glo)%npoints, 1, MPI_INT, domglo_rank(remote_ind_glo), & |
---|
296 | domloc_glo_ind(ind_loc)+ndomain_glo*remote_ind_glo, comm_icosa, reqs(ind_loc,remote_ind_glo,2), ierr ) |
---|
297 | end do |
---|
298 | end do |
---|
299 | call waitall(reqs,ndomain*ndomain_glo*2) |
---|
300 | |
---|
301 | ! Exchange points_send -> points_recv |
---|
302 | do ind_loc = 1, ndomain |
---|
303 | do remote_ind_glo = 1, ndomain_glo |
---|
304 | call MPI_Isend( points_send(ind_loc,remote_ind_glo)%i, points_send(ind_loc,remote_ind_glo)%npoints, MPI_INT, & |
---|
305 | domglo_rank(remote_ind_glo), remote_ind_glo+ndomain_glo*domloc_glo_ind(ind_loc), & |
---|
306 | comm_icosa, reqs(ind_loc,remote_ind_glo,1), ierr ) |
---|
307 | call MPI_Isend( points_send(ind_loc,remote_ind_glo)%j, points_send(ind_loc,remote_ind_glo)%npoints, MPI_INT, & |
---|
308 | domglo_rank(remote_ind_glo), remote_ind_glo+ndomain_glo*domloc_glo_ind(ind_loc), & |
---|
309 | comm_icosa, reqs2(ind_loc,remote_ind_glo,1), ierr ) |
---|
310 | |
---|
311 | points => req(ind_loc)%points_HtoB(remote_ind_glo) ! Points to recieve (HtoB) |
---|
312 | allocate( points%i( points%npoints ) ) |
---|
313 | allocate( points%j( points%npoints ) ) |
---|
314 | call MPI_Irecv( points%i, points%npoints, MPI_INT, domglo_rank(remote_ind_glo), & |
---|
315 | domloc_glo_ind(ind_loc)+ndomain_glo*remote_ind_glo, comm_icosa, reqs(ind_loc,remote_ind_glo,2), ierr) |
---|
316 | call MPI_Irecv( points%j, points%npoints, MPI_INT, domglo_rank(remote_ind_glo), & |
---|
317 | domloc_glo_ind(ind_loc)+ndomain_glo*remote_ind_glo, comm_icosa, reqs2(ind_loc,remote_ind_glo,2), ierr) |
---|
318 | |
---|
319 | if(req(1)%field_type == field_U) then |
---|
320 | call MPI_Isend( points_send(ind_loc,remote_ind_glo)%edge, points_send(ind_loc,remote_ind_glo)%npoints, MPI_INT, & |
---|
321 | domglo_rank(remote_ind_glo), remote_ind_glo+ndomain_glo*domloc_glo_ind(ind_loc), & |
---|
322 | comm_icosa, reqs3(ind_loc,remote_ind_glo,1), ierr ) |
---|
323 | allocate( points%edge( points%npoints ) ) |
---|
324 | call MPI_Irecv( points%edge, points%npoints, MPI_INT, & |
---|
325 | domglo_rank(remote_ind_glo), domloc_glo_ind(ind_loc)+ndomain_glo*remote_ind_glo, & |
---|
326 | comm_icosa, reqs3(ind_loc,remote_ind_glo,2), ierr ) |
---|
327 | endif |
---|
328 | end do |
---|
329 | end do |
---|
330 | call waitall(reqs,ndomain*ndomain_glo*2) |
---|
331 | call waitall(reqs2,ndomain*ndomain_glo*2) |
---|
332 | if(req(1)%field_type == field_U) call waitall(reqs3,ndomain*ndomain_glo*2) |
---|
333 | end subroutine |
---|
334 | |
---|
335 | ! MPI_Waitall, but with a N-Dim array of requests |
---|
336 | subroutine waitall(reqs, size) |
---|
337 | use mpi_mod |
---|
338 | integer :: size, reqs(size), ierr |
---|
339 | call MPI_Waitall( size, reqs, MPI_STATUSES_IGNORE, ierr) |
---|
340 | end subroutine |
---|
341 | end subroutine |
---|
342 | end module |
---|