1 | MODULE transfert_mpi_mod |
---|
2 | USE genmod |
---|
3 | |
---|
4 | TYPE array |
---|
5 | INTEGER,POINTER :: value(:) |
---|
6 | INTEGER,POINTER :: sign(:) |
---|
7 | INTEGER :: domain |
---|
8 | INTEGER :: rank |
---|
9 | INTEGER :: size |
---|
10 | INTEGER,POINTER :: buffer(:) |
---|
11 | REAL,POINTER :: buffer_r2(:) |
---|
12 | REAL,POINTER :: buffer_r3(:,:) |
---|
13 | REAL,POINTER :: buffer_r4(:,:,:) |
---|
14 | END TYPE array |
---|
15 | |
---|
16 | TYPE t_request |
---|
17 | INTEGER :: type_field |
---|
18 | INTEGER :: max_size |
---|
19 | INTEGER :: size |
---|
20 | LOGICAL :: vector |
---|
21 | INTEGER,POINTER :: src_domain(:) |
---|
22 | INTEGER,POINTER :: src_i(:) |
---|
23 | INTEGER,POINTER :: src_j(:) |
---|
24 | INTEGER,POINTER :: src_ind(:) |
---|
25 | INTEGER,POINTER :: target_ind(:) |
---|
26 | INTEGER,POINTER :: target_i(:) |
---|
27 | INTEGER,POINTER :: target_j(:) |
---|
28 | INTEGER,POINTER :: target_sign(:) |
---|
29 | INTEGER :: nrecv |
---|
30 | TYPE(ARRAY),POINTER :: recv(:) |
---|
31 | INTEGER :: nsend |
---|
32 | TYPE(ARRAY),POINTER :: send(:) |
---|
33 | END TYPE t_request |
---|
34 | |
---|
35 | TYPE(t_request),POINTER :: req_i1(:) |
---|
36 | TYPE(t_request),POINTER :: req_e1_scal(:) |
---|
37 | TYPE(t_request),POINTER :: req_e1_vect(:) |
---|
38 | |
---|
39 | TYPE(t_request),POINTER :: req_i0(:) |
---|
40 | TYPE(t_request),POINTER :: req_e0_scal(:) |
---|
41 | TYPE(t_request),POINTER :: req_e0_vect(:) |
---|
42 | |
---|
43 | |
---|
44 | CONTAINS |
---|
45 | |
---|
46 | SUBROUTINE init_transfert |
---|
47 | USE domain_mod |
---|
48 | USE dimensions |
---|
49 | USE field_mod |
---|
50 | USE metric |
---|
51 | USE mpipara |
---|
52 | IMPLICIT NONE |
---|
53 | INTEGER :: ind,i,j |
---|
54 | |
---|
55 | CALL create_request(field_t,req_i1) |
---|
56 | |
---|
57 | DO ind=1,ndomain |
---|
58 | CALL swap_dimensions(ind) |
---|
59 | DO i=ii_begin,ii_end+1 |
---|
60 | CALL request_add_point(ind,i,jj_begin-1,req_i1) |
---|
61 | ENDDO |
---|
62 | |
---|
63 | DO j=jj_begin,jj_end |
---|
64 | CALL request_add_point(ind,ii_end+1,j,req_i1) |
---|
65 | ENDDO |
---|
66 | DO i=ii_begin,ii_end |
---|
67 | CALL request_add_point(ind,i,jj_end+1,req_i1) |
---|
68 | ENDDO |
---|
69 | |
---|
70 | DO j=jj_begin,jj_end+1 |
---|
71 | CALL request_add_point(ind,ii_begin-1,j,req_i1) |
---|
72 | ENDDO |
---|
73 | |
---|
74 | DO i=ii_begin,ii_end |
---|
75 | ! CALL request_add_point(ind,i,jj_begin,req_i1) |
---|
76 | ENDDO |
---|
77 | |
---|
78 | DO j=jj_begin,jj_end |
---|
79 | ! CALL request_add_point(ind,ii_end,j,req_i1) |
---|
80 | ENDDO |
---|
81 | |
---|
82 | DO i=ii_begin,ii_end |
---|
83 | ! CALL request_add_point(ind,i,jj_end,req_i1) |
---|
84 | ENDDO |
---|
85 | |
---|
86 | DO j=jj_begin,jj_end |
---|
87 | ! CALL request_add_point(ind,ii_begin,j,req_i1) |
---|
88 | ENDDO |
---|
89 | |
---|
90 | ENDDO |
---|
91 | |
---|
92 | CALL finalize_request(req_i1) |
---|
93 | |
---|
94 | |
---|
95 | CALL create_request(field_t,req_i0) |
---|
96 | |
---|
97 | DO ind=1,ndomain |
---|
98 | CALL swap_dimensions(ind) |
---|
99 | |
---|
100 | DO i=ii_begin,ii_end |
---|
101 | CALL request_add_point(ind,i,jj_begin,req_i0) |
---|
102 | ENDDO |
---|
103 | |
---|
104 | DO j=jj_begin,jj_end |
---|
105 | CALL request_add_point(ind,ii_end,j,req_i0) |
---|
106 | ENDDO |
---|
107 | |
---|
108 | DO i=ii_begin,ii_end |
---|
109 | CALL request_add_point(ind,i,jj_end,req_i0) |
---|
110 | ENDDO |
---|
111 | |
---|
112 | DO j=jj_begin,jj_end |
---|
113 | CALL request_add_point(ind,ii_begin,j,req_i0) |
---|
114 | ENDDO |
---|
115 | |
---|
116 | ENDDO |
---|
117 | |
---|
118 | CALL finalize_request(req_i0) |
---|
119 | |
---|
120 | |
---|
121 | CALL create_request(field_u,req_e1_scal) |
---|
122 | DO ind=1,ndomain |
---|
123 | CALL swap_dimensions(ind) |
---|
124 | DO i=ii_begin,ii_end |
---|
125 | CALL request_add_point(ind,i,jj_begin-1,req_e1_scal,rup) |
---|
126 | CALL request_add_point(ind,i+1,jj_begin-1,req_e1_scal,lup) |
---|
127 | ENDDO |
---|
128 | |
---|
129 | DO j=jj_begin,jj_end |
---|
130 | CALL request_add_point(ind,ii_end+1,j,req_e1_scal,left) |
---|
131 | CALL request_add_point(ind,ii_end+1,j-1,req_e1_scal,lup) |
---|
132 | ENDDO |
---|
133 | |
---|
134 | DO i=ii_begin,ii_end |
---|
135 | CALL request_add_point(ind,i,jj_end+1,req_e1_scal,ldown) |
---|
136 | CALL request_add_point(ind,i-1,jj_end+1,req_e1_scal,rdown) |
---|
137 | ENDDO |
---|
138 | |
---|
139 | DO j=jj_begin,jj_end |
---|
140 | CALL request_add_point(ind,ii_begin-1,j,req_e1_scal,right) |
---|
141 | CALL request_add_point(ind,ii_begin-1,j+1,req_e1_scal,rdown) |
---|
142 | ENDDO |
---|
143 | |
---|
144 | DO i=ii_begin+1,ii_end-1 |
---|
145 | ! CALL request_add_point(ind,i,jj_begin,req_e1_scal,right) |
---|
146 | ! CALL request_add_point(ind,i,jj_end,req_e1_scal,right) |
---|
147 | ENDDO |
---|
148 | |
---|
149 | DO j=jj_begin+1,jj_end-1 |
---|
150 | ! CALL request_add_point(ind,ii_begin,j,req_e1_scal,rup) |
---|
151 | ! CALL request_add_point(ind,ii_end,j,req_e1_scal,rup) |
---|
152 | ENDDO |
---|
153 | |
---|
154 | ! CALL request_add_point(ind,ii_begin+1,jj_begin,req_e1_scal,left) |
---|
155 | ! CALL request_add_point(ind,ii_begin,jj_begin+1,req_e1_scal,ldown) |
---|
156 | ! CALL request_add_point(ind,ii_begin+1,jj_end,req_e1_scal,left) |
---|
157 | ! CALL request_add_point(ind,ii_end,jj_begin+1,req_e1_scal,ldown) |
---|
158 | |
---|
159 | ENDDO |
---|
160 | |
---|
161 | CALL finalize_request(req_e1_scal) |
---|
162 | |
---|
163 | |
---|
164 | CALL create_request(field_u,req_e0_scal) |
---|
165 | DO ind=1,ndomain |
---|
166 | CALL swap_dimensions(ind) |
---|
167 | |
---|
168 | |
---|
169 | DO i=ii_begin+1,ii_end-1 |
---|
170 | CALL request_add_point(ind,i,jj_begin,req_e0_scal,right) |
---|
171 | CALL request_add_point(ind,i,jj_end,req_e0_scal,right) |
---|
172 | ENDDO |
---|
173 | |
---|
174 | DO j=jj_begin+1,jj_end-1 |
---|
175 | CALL request_add_point(ind,ii_begin,j,req_e0_scal,rup) |
---|
176 | CALL request_add_point(ind,ii_end,j,req_e0_scal,rup) |
---|
177 | ENDDO |
---|
178 | |
---|
179 | CALL request_add_point(ind,ii_begin+1,jj_begin,req_e0_scal,left) |
---|
180 | CALL request_add_point(ind,ii_begin,jj_begin+1,req_e0_scal,ldown) |
---|
181 | CALL request_add_point(ind,ii_begin+1,jj_end,req_e0_scal,left) |
---|
182 | CALL request_add_point(ind,ii_end,jj_begin+1,req_e0_scal,ldown) |
---|
183 | |
---|
184 | ENDDO |
---|
185 | |
---|
186 | CALL finalize_request(req_e0_scal) |
---|
187 | |
---|
188 | |
---|
189 | |
---|
190 | CALL create_request(field_u,req_e1_vect,.TRUE.) |
---|
191 | DO ind=1,ndomain |
---|
192 | CALL swap_dimensions(ind) |
---|
193 | DO i=ii_begin,ii_end |
---|
194 | CALL request_add_point(ind,i,jj_begin-1,req_e1_vect,rup) |
---|
195 | CALL request_add_point(ind,i+1,jj_begin-1,req_e1_vect,lup) |
---|
196 | ENDDO |
---|
197 | |
---|
198 | DO j=jj_begin,jj_end |
---|
199 | CALL request_add_point(ind,ii_end+1,j,req_e1_vect,left) |
---|
200 | CALL request_add_point(ind,ii_end+1,j-1,req_e1_vect,lup) |
---|
201 | ENDDO |
---|
202 | |
---|
203 | DO i=ii_begin,ii_end |
---|
204 | CALL request_add_point(ind,i,jj_end+1,req_e1_vect,ldown) |
---|
205 | CALL request_add_point(ind,i-1,jj_end+1,req_e1_vect,rdown) |
---|
206 | ENDDO |
---|
207 | |
---|
208 | DO j=jj_begin,jj_end |
---|
209 | CALL request_add_point(ind,ii_begin-1,j,req_e1_vect,right) |
---|
210 | CALL request_add_point(ind,ii_begin-1,j+1,req_e1_vect,rdown) |
---|
211 | ENDDO |
---|
212 | |
---|
213 | DO i=ii_begin+1,ii_end-1 |
---|
214 | ! CALL request_add_point(ind,i,jj_begin,req_e1_vect,right) |
---|
215 | ! CALL request_add_point(ind,i,jj_end,req_e1_vect,right) |
---|
216 | ENDDO |
---|
217 | |
---|
218 | DO j=jj_begin+1,jj_end-1 |
---|
219 | ! CALL request_add_point(ind,ii_begin,j,req_e1_vect,rup) |
---|
220 | ! CALL request_add_point(ind,ii_end,j,req_e1_vect,rup) |
---|
221 | ENDDO |
---|
222 | |
---|
223 | ! CALL request_add_point(ind,ii_begin+1,jj_begin,req_e1_vect,left) |
---|
224 | ! CALL request_add_point(ind,ii_begin,jj_begin+1,req_e1_vect,ldown) |
---|
225 | ! CALL request_add_point(ind,ii_begin+1,jj_end,req_e1_vect,left) |
---|
226 | ! CALL request_add_point(ind,ii_end,jj_begin+1,req_e1_vect,ldown) |
---|
227 | |
---|
228 | ENDDO |
---|
229 | |
---|
230 | CALL finalize_request(req_e1_vect) |
---|
231 | |
---|
232 | |
---|
233 | CALL create_request(field_u,req_e0_vect,.TRUE.) |
---|
234 | DO ind=1,ndomain |
---|
235 | CALL swap_dimensions(ind) |
---|
236 | |
---|
237 | DO i=ii_begin+1,ii_end-1 |
---|
238 | CALL request_add_point(ind,i,jj_begin,req_e0_vect,right) |
---|
239 | CALL request_add_point(ind,i,jj_end,req_e0_vect,right) |
---|
240 | ENDDO |
---|
241 | |
---|
242 | DO j=jj_begin+1,jj_end-1 |
---|
243 | CALL request_add_point(ind,ii_begin,j,req_e0_vect,rup) |
---|
244 | CALL request_add_point(ind,ii_end,j,req_e0_vect,rup) |
---|
245 | ENDDO |
---|
246 | |
---|
247 | CALL request_add_point(ind,ii_begin+1,jj_begin,req_e0_vect,left) |
---|
248 | CALL request_add_point(ind,ii_begin,jj_begin+1,req_e0_vect,ldown) |
---|
249 | CALL request_add_point(ind,ii_begin+1,jj_end,req_e0_vect,left) |
---|
250 | CALL request_add_point(ind,ii_end,jj_begin+1,req_e0_vect,ldown) |
---|
251 | |
---|
252 | ENDDO |
---|
253 | |
---|
254 | CALL finalize_request(req_e0_vect) |
---|
255 | |
---|
256 | |
---|
257 | END SUBROUTINE init_transfert |
---|
258 | |
---|
259 | SUBROUTINE create_request(type_field,request,vector) |
---|
260 | USE domain_mod |
---|
261 | USE field_mod |
---|
262 | IMPLICIT NONE |
---|
263 | INTEGER :: type_field |
---|
264 | TYPE(t_request),POINTER :: request(:) |
---|
265 | LOGICAL,OPTIONAL :: vector |
---|
266 | |
---|
267 | TYPE(t_request),POINTER :: req |
---|
268 | TYPE(t_domain),POINTER :: d |
---|
269 | INTEGER :: ind |
---|
270 | INTEGER :: max_size |
---|
271 | |
---|
272 | ALLOCATE(request(ndomain)) |
---|
273 | |
---|
274 | DO ind=1,ndomain |
---|
275 | req=>request(ind) |
---|
276 | d=>domain(ind) |
---|
277 | IF (type_field==field_t) THEN |
---|
278 | Max_size=2*(d%iim+2)+2*(d%jjm+2) |
---|
279 | ELSE IF (type_field==field_u) THEN |
---|
280 | Max_size=3*(2*(d%iim+2)+2*(d%jjm+2)) |
---|
281 | ELSE IF (type_field==field_z) THEN |
---|
282 | Max_size=2*(2*(d%iim+2)+2*(d%jjm+2)) |
---|
283 | ENDIF |
---|
284 | |
---|
285 | req%type_field=type_field |
---|
286 | req%max_size=max_size*2 |
---|
287 | req%size=0 |
---|
288 | req%vector=.FALSE. |
---|
289 | IF (PRESENT(vector)) req%vector=vector |
---|
290 | ALLOCATE(req%src_domain(req%max_size)) |
---|
291 | ALLOCATE(req%src_ind(req%max_size)) |
---|
292 | ALLOCATE(req%target_ind(req%max_size)) |
---|
293 | ALLOCATE(req%src_i(req%max_size)) |
---|
294 | ALLOCATE(req%src_j(req%max_size)) |
---|
295 | ALLOCATE(req%target_i(req%max_size)) |
---|
296 | ALLOCATE(req%target_j(req%max_size)) |
---|
297 | ALLOCATE(req%target_sign(req%max_size)) |
---|
298 | ENDDO |
---|
299 | |
---|
300 | END SUBROUTINE create_request |
---|
301 | |
---|
302 | SUBROUTINE reallocate_request(req) |
---|
303 | IMPLICIT NONE |
---|
304 | TYPE(t_request),POINTER :: req |
---|
305 | |
---|
306 | INTEGER,POINTER :: src_domain(:) |
---|
307 | INTEGER,POINTER :: src_ind(:) |
---|
308 | INTEGER,POINTER :: target_ind(:) |
---|
309 | INTEGER,POINTER :: src_i(:) |
---|
310 | INTEGER,POINTER :: src_j(:) |
---|
311 | INTEGER,POINTER :: target_i(:) |
---|
312 | INTEGER,POINTER :: target_j(:) |
---|
313 | INTEGER,POINTER :: target_sign(:) |
---|
314 | |
---|
315 | PRINT *,"REALLOCATE_REQUEST" |
---|
316 | src_domain=>req%src_domain |
---|
317 | src_ind=>req%src_ind |
---|
318 | target_ind=>req%target_ind |
---|
319 | src_i=>req%src_i |
---|
320 | src_j=>req%src_j |
---|
321 | target_i=>req%target_i |
---|
322 | target_j=>req%target_j |
---|
323 | target_sign=>req%target_sign |
---|
324 | ! req%max_size=req%max_size*2 |
---|
325 | ALLOCATE(req%src_domain(req%max_size*2)) |
---|
326 | ALLOCATE(req%src_ind(req%max_size*2)) |
---|
327 | ALLOCATE(req%target_ind(req%max_size*2)) |
---|
328 | ALLOCATE(req%src_i(req%max_size*2)) |
---|
329 | ALLOCATE(req%src_j(req%max_size*2)) |
---|
330 | ALLOCATE(req%target_i(req%max_size*2)) |
---|
331 | ALLOCATE(req%target_j(req%max_size*2)) |
---|
332 | ALLOCATE(req%target_sign(req%max_size*2)) |
---|
333 | |
---|
334 | req%src_domain(1:req%max_size)=src_domain(:) |
---|
335 | req%src_ind(1:req%max_size)=src_ind(:) |
---|
336 | req%target_ind(1:req%max_size)=target_ind(:) |
---|
337 | req%src_i(1:req%max_size)=src_i(:) |
---|
338 | req%src_j(1:req%max_size)=src_j(:) |
---|
339 | req%target_i(1:req%max_size)=target_i(:) |
---|
340 | req%target_j(1:req%max_size)=target_j(:) |
---|
341 | req%target_sign(1:req%max_size)=target_sign(:) |
---|
342 | |
---|
343 | req%max_size=req%max_size*2 |
---|
344 | |
---|
345 | DEALLOCATE(src_domain) |
---|
346 | DEALLOCATE(src_ind) |
---|
347 | DEALLOCATE(target_ind) |
---|
348 | DEALLOCATE(src_i) |
---|
349 | DEALLOCATE(src_j) |
---|
350 | DEALLOCATE(target_i) |
---|
351 | DEALLOCATE(target_j) |
---|
352 | DEALLOCATE(target_sign) |
---|
353 | |
---|
354 | END SUBROUTINE reallocate_request |
---|
355 | |
---|
356 | |
---|
357 | SUBROUTINE request_add_point(ind,i,j,request,pos) |
---|
358 | USE domain_mod |
---|
359 | USE field_mod |
---|
360 | IMPLICIT NONE |
---|
361 | INTEGER,INTENT(IN) :: ind |
---|
362 | INTEGER,INTENT(IN) :: i |
---|
363 | INTEGER,INTENT(IN) :: j |
---|
364 | TYPE(t_request),POINTER :: request(:) |
---|
365 | INTEGER,INTENT(IN),OPTIONAL :: pos |
---|
366 | |
---|
367 | INTEGER :: src_domain |
---|
368 | INTEGER :: src_iim,src_i,src_j,src_n,src_pos,src_delta |
---|
369 | TYPE(t_request),POINTER :: req |
---|
370 | TYPE(t_domain),POINTER :: d |
---|
371 | |
---|
372 | req=>request(ind) |
---|
373 | d=>domain(ind) |
---|
374 | |
---|
375 | IF (req%max_size==req%size) CALL reallocate_request(req) |
---|
376 | req%size=req%size+1 |
---|
377 | IF (req%type_field==field_t) THEN |
---|
378 | src_domain=domain(ind)%assign_domain(i,j) |
---|
379 | src_iim=domain_glo(src_domain)%iim |
---|
380 | src_i=domain(ind)%assign_i(i,j) |
---|
381 | src_j=domain(ind)%assign_j(i,j) |
---|
382 | |
---|
383 | req%target_ind(req%size)=(j-1)*d%iim+i |
---|
384 | req%target_sign(req%size)=1 |
---|
385 | req%src_domain(req%size)=src_domain |
---|
386 | req%src_ind(req%size)=(src_j-1)*src_iim+src_i |
---|
387 | ELSE IF (req%type_field==field_u) THEN |
---|
388 | IF (.NOT. PRESENT(pos)) STOP 'argument request_add_point non conforme' |
---|
389 | |
---|
390 | src_domain=domain(ind)%edge_assign_domain(pos-1,i,j) |
---|
391 | src_iim=domain_glo(src_domain)%iim |
---|
392 | src_i=domain(ind)%edge_assign_i(pos-1,i,j) |
---|
393 | src_j=domain(ind)%edge_assign_j(pos-1,i,j) |
---|
394 | src_n=(src_j-1)*src_iim+src_i |
---|
395 | src_delta=domain(ind)%delta(i,j) |
---|
396 | src_pos=domain(ind)%edge_assign_pos(pos-1,i,j)+1 |
---|
397 | |
---|
398 | req%target_ind(req%size)=(j-1)*d%iim+i+d%u_pos(pos) |
---|
399 | |
---|
400 | req%target_sign(req%size)= 1 |
---|
401 | IF (req%vector) req%target_sign(req%size)= domain(ind)%edge_assign_sign(pos-1,i,j) |
---|
402 | |
---|
403 | req%src_domain(req%size)=src_domain |
---|
404 | req%src_ind(req%size)=src_n+domain_glo(src_domain)%u_pos(src_pos) |
---|
405 | |
---|
406 | ELSE IF (req%type_field==field_z) THEN |
---|
407 | IF (.NOT. PRESENT(pos)) STOP 'argument request_add_point non conforme' |
---|
408 | |
---|
409 | src_domain=domain(ind)%assign_domain(i,j) |
---|
410 | src_iim=domain_glo(src_domain)%iim |
---|
411 | src_i=domain(ind)%assign_i(i,j) |
---|
412 | src_j=domain(ind)%assign_j(i,j) |
---|
413 | src_n=(src_j-1)*src_iim+src_i |
---|
414 | src_delta=domain(ind)%delta(i,j) |
---|
415 | |
---|
416 | src_pos=MOD(pos-1+src_delta+6,6)+1 |
---|
417 | |
---|
418 | req%target_ind(req%size)=(j-1)*d%iim+i+d%z_pos(pos) |
---|
419 | req%target_sign(req%size)=1 |
---|
420 | req%src_domain(req%size)=src_domain |
---|
421 | req%src_ind(req%size)=src_n+domain_glo(src_domain)%z_pos(src_pos) |
---|
422 | ENDIF |
---|
423 | END SUBROUTINE request_add_point |
---|
424 | |
---|
425 | |
---|
426 | SUBROUTINE Finalize_request(request) |
---|
427 | USE mpipara |
---|
428 | USE domain_mod |
---|
429 | USE mpi_mod |
---|
430 | IMPLICIT NONE |
---|
431 | TYPE(t_request),POINTER :: request(:) |
---|
432 | TYPE(t_request),POINTER :: req |
---|
433 | INTEGER :: nb_domain_recv(0:mpi_size-1) |
---|
434 | INTEGER :: nb_domain_send(0:mpi_size-1) |
---|
435 | INTEGER :: nb_data_domain_recv(ndomain_glo) |
---|
436 | INTEGER :: list_domain_recv(ndomain_glo) |
---|
437 | INTEGER,ALLOCATABLE :: list_domain_send(:) |
---|
438 | INTEGER :: list_domain(ndomain) |
---|
439 | |
---|
440 | INTEGER :: rank,i,j |
---|
441 | INTEGER :: size,ind_glo,ind_loc |
---|
442 | INTEGER :: isend, irecv, ireq, nreq |
---|
443 | INTEGER, ALLOCATABLE :: mpi_req(:) |
---|
444 | INTEGER, ALLOCATABLE :: status(:,:) |
---|
445 | |
---|
446 | IF (.NOT. using_mpi) RETURN |
---|
447 | |
---|
448 | DO ind_loc=1,ndomain |
---|
449 | req=>request(ind_loc) |
---|
450 | |
---|
451 | nb_data_domain_recv(:) = 0 |
---|
452 | nb_domain_recv(:) = 0 |
---|
453 | |
---|
454 | DO i=1,req%size |
---|
455 | ind_glo=req%src_domain(i) |
---|
456 | nb_data_domain_recv(ind_glo)=nb_data_domain_recv(ind_glo)+1 |
---|
457 | ENDDO |
---|
458 | |
---|
459 | DO ind_glo=1,ndomain_glo |
---|
460 | IF ( nb_data_domain_recv(ind_glo) > 0 ) nb_domain_recv(domglo_rank(ind_glo))=nb_domain_recv(domglo_rank(ind_glo))+1 |
---|
461 | ENDDO |
---|
462 | |
---|
463 | req%nrecv=sum(nb_domain_recv(:)) |
---|
464 | ALLOCATE(req%recv(req%nrecv)) |
---|
465 | |
---|
466 | irecv=0 |
---|
467 | DO ind_glo=1,ndomain_glo |
---|
468 | IF (nb_data_domain_recv(ind_glo)>0) THEN |
---|
469 | irecv=irecv+1 |
---|
470 | list_domain_recv(ind_glo)=irecv |
---|
471 | req%recv(irecv)%rank=domglo_rank(ind_glo) |
---|
472 | req%recv(irecv)%size=nb_data_domain_recv(ind_glo) |
---|
473 | req%recv(irecv)%domain=domglo_loc_ind(ind_glo) |
---|
474 | ALLOCATE(req%recv(irecv)%value(req%recv(irecv)%size)) |
---|
475 | ALLOCATE(req%recv(irecv)%sign(req%recv(irecv)%size)) |
---|
476 | ALLOCATE(req%recv(irecv)%buffer(req%recv(irecv)%size)) |
---|
477 | ENDIF |
---|
478 | ENDDO |
---|
479 | |
---|
480 | req%recv(:)%size=0 |
---|
481 | irecv=0 |
---|
482 | DO i=1,req%size |
---|
483 | irecv=list_domain_recv(req%src_domain(i)) |
---|
484 | req%recv(irecv)%size=req%recv(irecv)%size+1 |
---|
485 | size=req%recv(irecv)%size |
---|
486 | req%recv(irecv)%value(size)=req%src_ind(i) |
---|
487 | req%recv(irecv)%buffer(size)=req%target_ind(i) |
---|
488 | req%recv(irecv)%sign(size)=req%target_sign(i) |
---|
489 | ENDDO |
---|
490 | ENDDO |
---|
491 | |
---|
492 | nb_domain_recv(:) = 0 |
---|
493 | DO ind_loc=1,ndomain |
---|
494 | req=>request(ind_loc) |
---|
495 | |
---|
496 | DO irecv=1,req%nrecv |
---|
497 | rank=req%recv(irecv)%rank |
---|
498 | nb_domain_recv(rank)=nb_domain_recv(rank)+1 |
---|
499 | ENDDO |
---|
500 | ENDDO |
---|
501 | |
---|
502 | CALL MPI_Alltoall(nb_domain_recv,1,MPI_INTEGER,nb_domain_send,1,MPI_INTEGER,comm_icosa,ierr) |
---|
503 | |
---|
504 | |
---|
505 | ALLOCATE(list_domain_send(sum(nb_domain_send))) |
---|
506 | |
---|
507 | nreq=sum(nb_domain_recv(:))+sum(nb_domain_send(:)) |
---|
508 | ALLOCATE(mpi_req(nreq)) |
---|
509 | ALLOCATE(status(MPI_STATUS_SIZE,nreq)) |
---|
510 | |
---|
511 | ireq=0 |
---|
512 | DO ind_loc=1,ndomain |
---|
513 | req=>request(ind_loc) |
---|
514 | DO irecv=1,req%nrecv |
---|
515 | ireq=ireq+1 |
---|
516 | CALL MPI_ISEND(req%recv(irecv)%domain,1,MPI_INTEGER,req%recv(irecv)%rank,0,comm_icosa, mpi_req(ireq),ierr) |
---|
517 | ENDDO |
---|
518 | ENDDO |
---|
519 | |
---|
520 | j=0 |
---|
521 | DO rank=0,mpi_size-1 |
---|
522 | DO i=1,nb_domain_send(rank) |
---|
523 | j=j+1 |
---|
524 | ireq=ireq+1 |
---|
525 | CALL MPI_IRECV(list_domain_send(j),1,MPI_INTEGER,rank,0,comm_icosa, mpi_req(ireq),ierr) |
---|
526 | ENDDO |
---|
527 | ENDDO |
---|
528 | |
---|
529 | CALL MPI_WAITALL(nreq,mpi_req,status,ierr) |
---|
530 | |
---|
531 | list_domain(:)=0 |
---|
532 | DO i=1,sum(nb_domain_send) |
---|
533 | ind_loc=list_domain_send(i) |
---|
534 | list_domain(ind_loc)=list_domain(ind_loc)+1 |
---|
535 | ENDDO |
---|
536 | |
---|
537 | DO ind_loc=1,ndomain |
---|
538 | req=>request(ind_loc) |
---|
539 | req%nsend=list_domain(ind_loc) |
---|
540 | ALLOCATE(req%send(req%nsend)) |
---|
541 | ENDDO |
---|
542 | |
---|
543 | ireq=0 |
---|
544 | DO ind_loc=1,ndomain |
---|
545 | req=>request(ind_loc) |
---|
546 | |
---|
547 | DO irecv=1,req%nrecv |
---|
548 | ireq=ireq+1 |
---|
549 | CALL MPI_ISEND(mpi_rank,1,MPI_INTEGER,req%recv(irecv)%rank,req%recv(irecv)%domain,comm_icosa, mpi_req(ireq),ierr) |
---|
550 | ENDDO |
---|
551 | |
---|
552 | DO isend=1,req%nsend |
---|
553 | ireq=ireq+1 |
---|
554 | CALL MPI_IRECV(req%send(isend)%rank,1,MPI_INTEGER,MPI_ANY_SOURCE,ind_loc,comm_icosa, mpi_req(ireq),ierr) |
---|
555 | ENDDO |
---|
556 | ENDDO |
---|
557 | |
---|
558 | CALL MPI_WAITALL(nreq,mpi_req,status,ierr) |
---|
559 | CALL MPI_BARRIER(comm_icosa,ierr) |
---|
560 | |
---|
561 | ireq=0 |
---|
562 | DO ind_loc=1,ndomain |
---|
563 | req=>request(ind_loc) |
---|
564 | |
---|
565 | DO irecv=1,req%nrecv |
---|
566 | ireq=ireq+1 |
---|
567 | CALL MPI_ISEND(req%recv(irecv)%size,1,MPI_INTEGER,req%recv(irecv)%rank,req%recv(irecv)%domain,comm_icosa, mpi_req(ireq),ierr) |
---|
568 | ENDDO |
---|
569 | |
---|
570 | DO isend=1,req%nsend |
---|
571 | ireq=ireq+1 |
---|
572 | CALL MPI_IRECV(req%send(isend)%size,1,MPI_INTEGER,req%send(isend)%rank,ind_loc,comm_icosa, mpi_req(ireq),ierr) |
---|
573 | ENDDO |
---|
574 | ENDDO |
---|
575 | |
---|
576 | CALL MPI_WAITALL(nreq,mpi_req,status,ierr) |
---|
577 | |
---|
578 | ireq=0 |
---|
579 | DO ind_loc=1,ndomain |
---|
580 | req=>request(ind_loc) |
---|
581 | |
---|
582 | DO irecv=1,req%nrecv |
---|
583 | ireq=ireq+1 |
---|
584 | CALL MPI_ISEND(req%recv(irecv)%value,req%recv(irecv)%size,MPI_INTEGER,& |
---|
585 | req%recv(irecv)%rank,req%recv(irecv)%domain,comm_icosa, mpi_req(ireq),ierr) |
---|
586 | ENDDO |
---|
587 | |
---|
588 | DO isend=1,req%nsend |
---|
589 | ireq=ireq+1 |
---|
590 | ALLOCATE(req%send(isend)%value(req%send(isend)%size)) |
---|
591 | CALL MPI_IRECV(req%send(isend)%value,req%send(isend)%size,MPI_INTEGER,& |
---|
592 | req%send(isend)%rank,ind_loc,comm_icosa, mpi_req(ireq),ierr) |
---|
593 | ENDDO |
---|
594 | ENDDO |
---|
595 | |
---|
596 | CALL MPI_WAITALL(nreq,mpi_req,status,ierr) |
---|
597 | |
---|
598 | DO ind_loc=1,ndomain |
---|
599 | req=>request(ind_loc) |
---|
600 | |
---|
601 | DO irecv=1,req%nrecv |
---|
602 | req%recv(irecv)%value(:)=req%recv(irecv)%buffer(:) |
---|
603 | req%recv(irecv)%sign(:) =req%recv(irecv)%sign(:) |
---|
604 | DEALLOCATE(req%recv(irecv)%buffer) |
---|
605 | ENDDO |
---|
606 | ENDDO |
---|
607 | |
---|
608 | |
---|
609 | END SUBROUTINE Finalize_request |
---|
610 | |
---|
611 | |
---|
612 | SUBROUTINE transfert_request_mpi(field,request) |
---|
613 | USE field_mod |
---|
614 | USE domain_mod |
---|
615 | USE mpi_mod |
---|
616 | USE mpipara |
---|
617 | USE trace |
---|
618 | IMPLICIT NONE |
---|
619 | TYPE(t_field),POINTER :: field(:) |
---|
620 | TYPE(t_request),POINTER :: request(:) |
---|
621 | REAL(rstd),POINTER :: rval2d(:) |
---|
622 | REAL(rstd),POINTER :: rval3d(:,:) |
---|
623 | REAL(rstd),POINTER :: rval4d(:,:,:) |
---|
624 | REAL(rstd),POINTER :: buffer_r2(:) |
---|
625 | REAL(rstd),POINTER :: buffer_r3(:,:) |
---|
626 | REAL(rstd),POINTER :: buffer_r4(:,:,:) |
---|
627 | INTEGER,POINTER :: value(:) |
---|
628 | INTEGER,POINTER :: sgn(:) |
---|
629 | TYPE(ARRAY),POINTER :: recv,send |
---|
630 | TYPE(t_request),POINTER :: req |
---|
631 | INTEGER, ALLOCATABLE :: mpi_req(:) |
---|
632 | INTEGER, ALLOCATABLE :: status(:,:) |
---|
633 | INTEGER :: irecv,isend |
---|
634 | INTEGER :: ireq,nreq |
---|
635 | INTEGER :: ind,n |
---|
636 | INTEGER :: dim3,dim4 |
---|
637 | |
---|
638 | CALL trace_start("transfert_mpi") |
---|
639 | |
---|
640 | IF (field(1)%data_type==type_real) THEN |
---|
641 | IF (field(1)%ndim==2) THEN |
---|
642 | |
---|
643 | nreq=sum(request(:)%nsend)+sum(request(:)%nrecv) |
---|
644 | ALLOCATE(mpi_req(nreq)) |
---|
645 | ALLOCATE(status(MPI_STATUS_SIZE,nreq)) |
---|
646 | |
---|
647 | ireq=0 |
---|
648 | DO ind=1,ndomain |
---|
649 | rval2d=>field(ind)%rval2d |
---|
650 | |
---|
651 | req=>request(ind) |
---|
652 | DO isend=1,req%nsend |
---|
653 | send=>req%send(isend) |
---|
654 | |
---|
655 | ALLOCATE(send%buffer_r2(send%size)) |
---|
656 | buffer_r2=>send%buffer_r2 |
---|
657 | value=>send%value |
---|
658 | DO n=1,send%size |
---|
659 | buffer_r2(n)=rval2d(value(n)) |
---|
660 | ENDDO |
---|
661 | |
---|
662 | ireq=ireq+1 |
---|
663 | CALL MPI_ISEND(buffer_r2,send%size,MPI_REAL8,send%rank,ind,comm_icosa, mpi_req(ireq),ierr) |
---|
664 | ENDDO |
---|
665 | |
---|
666 | DO irecv=1,req%nrecv |
---|
667 | recv=>req%recv(irecv) |
---|
668 | ALLOCATE(recv%buffer_r2(recv%size)) |
---|
669 | |
---|
670 | ireq=ireq+1 |
---|
671 | CALL MPI_IRECV(recv%buffer_r2,recv%size,MPI_REAL8,recv%rank,recv%domain,comm_icosa, mpi_req(ireq),ierr) |
---|
672 | ENDDO |
---|
673 | |
---|
674 | ENDDO |
---|
675 | |
---|
676 | CALL MPI_WAITALL(nreq,mpi_req,status,ierr) |
---|
677 | ! DO ind=1,ndomain |
---|
678 | ! field(ind)%rval2d(:)=0. |
---|
679 | ! ENDDO |
---|
680 | |
---|
681 | DO ind=1,ndomain |
---|
682 | rval2d=>field(ind)%rval2d |
---|
683 | |
---|
684 | req=>request(ind) |
---|
685 | DO isend=1,req%nsend |
---|
686 | send=>req%send(isend) |
---|
687 | DEALLOCATE(send%buffer_r2) |
---|
688 | ENDDO |
---|
689 | |
---|
690 | DO irecv=1,req%nrecv |
---|
691 | recv=>req%recv(irecv) |
---|
692 | buffer_r2=>recv%buffer_r2 |
---|
693 | value=>recv%value |
---|
694 | sgn=>recv%sign |
---|
695 | DO n=1,recv%size |
---|
696 | rval2d(value(n))=buffer_r2(n)*sgn(n) |
---|
697 | ENDDO |
---|
698 | DEALLOCATE(recv%buffer_r2) |
---|
699 | ENDDO |
---|
700 | |
---|
701 | ENDDO |
---|
702 | |
---|
703 | |
---|
704 | ELSE IF (field(1)%ndim==3) THEN |
---|
705 | |
---|
706 | nreq=sum(request(:)%nsend)+sum(request(:)%nrecv) |
---|
707 | ALLOCATE(mpi_req(nreq)) |
---|
708 | ALLOCATE(status(MPI_STATUS_SIZE,nreq)) |
---|
709 | |
---|
710 | ireq=0 |
---|
711 | DO ind=1,ndomain |
---|
712 | dim3=size(field(ind)%rval3d,2) |
---|
713 | rval3d=>field(ind)%rval3d |
---|
714 | |
---|
715 | req=>request(ind) |
---|
716 | DO isend=1,req%nsend |
---|
717 | send=>req%send(isend) |
---|
718 | |
---|
719 | ALLOCATE(send%buffer_r3(send%size,dim3)) |
---|
720 | buffer_r3=>send%buffer_r3 |
---|
721 | value=>send%value |
---|
722 | DO n=1,send%size |
---|
723 | buffer_r3(n,:)=rval3d(value(n),:) |
---|
724 | ENDDO |
---|
725 | |
---|
726 | ireq=ireq+1 |
---|
727 | CALL MPI_ISEND(buffer_r3,send%size*dim3,MPI_REAL8,send%rank,ind,comm_icosa, mpi_req(ireq),ierr) |
---|
728 | ENDDO |
---|
729 | |
---|
730 | DO irecv=1,req%nrecv |
---|
731 | recv=>req%recv(irecv) |
---|
732 | ALLOCATE(recv%buffer_r3(recv%size,dim3)) |
---|
733 | |
---|
734 | ireq=ireq+1 |
---|
735 | CALL MPI_IRECV(recv%buffer_r3,recv%size*dim3,MPI_REAL8,recv%rank,recv%domain,comm_icosa, mpi_req(ireq),ierr) |
---|
736 | ENDDO |
---|
737 | |
---|
738 | ENDDO |
---|
739 | |
---|
740 | CALL MPI_WAITALL(nreq,mpi_req,status,ierr) |
---|
741 | ! DO ind=1,ndomain |
---|
742 | ! field(ind)%rval2d(:)=0. |
---|
743 | ! ENDDO |
---|
744 | |
---|
745 | DO ind=1,ndomain |
---|
746 | rval3d=>field(ind)%rval3d |
---|
747 | |
---|
748 | req=>request(ind) |
---|
749 | DO isend=1,req%nsend |
---|
750 | send=>req%send(isend) |
---|
751 | DEALLOCATE(send%buffer_r3) |
---|
752 | ENDDO |
---|
753 | |
---|
754 | DO irecv=1,req%nrecv |
---|
755 | recv=>req%recv(irecv) |
---|
756 | buffer_r3=>recv%buffer_r3 |
---|
757 | value=>recv%value |
---|
758 | sgn=>recv%sign |
---|
759 | DO n=1,recv%size |
---|
760 | rval3d(value(n),:)=buffer_r3(n,:)*sgn(n) |
---|
761 | ENDDO |
---|
762 | DEALLOCATE(recv%buffer_r3) |
---|
763 | ENDDO |
---|
764 | |
---|
765 | ENDDO |
---|
766 | |
---|
767 | ELSE IF (field(1)%ndim==4) THEN |
---|
768 | |
---|
769 | nreq=sum(request(:)%nsend)+sum(request(:)%nrecv) |
---|
770 | ALLOCATE(mpi_req(nreq)) |
---|
771 | ALLOCATE(status(MPI_STATUS_SIZE,nreq)) |
---|
772 | |
---|
773 | ireq=0 |
---|
774 | DO ind=1,ndomain |
---|
775 | dim3=size(field(ind)%rval4d,2) |
---|
776 | dim4=size(field(ind)%rval4d,3) |
---|
777 | rval4d=>field(ind)%rval4d |
---|
778 | |
---|
779 | req=>request(ind) |
---|
780 | DO isend=1,req%nsend |
---|
781 | send=>req%send(isend) |
---|
782 | |
---|
783 | ALLOCATE(send%buffer_r4(send%size,dim3,dim4)) |
---|
784 | buffer_r4=>send%buffer_r4 |
---|
785 | value=>send%value |
---|
786 | DO n=1,send%size |
---|
787 | buffer_r4(n,:,:)=rval4d(value(n),:,:) |
---|
788 | ENDDO |
---|
789 | |
---|
790 | ireq=ireq+1 |
---|
791 | CALL MPI_ISEND(buffer_r4,send%size*dim3*dim4,MPI_REAL8,send%rank,ind,comm_icosa, mpi_req(ireq),ierr) |
---|
792 | ENDDO |
---|
793 | |
---|
794 | DO irecv=1,req%nrecv |
---|
795 | recv=>req%recv(irecv) |
---|
796 | ALLOCATE(recv%buffer_r4(recv%size,dim3,dim4)) |
---|
797 | |
---|
798 | ireq=ireq+1 |
---|
799 | CALL MPI_IRECV(recv%buffer_r4,recv%size*dim3*dim4,MPI_REAL8,recv%rank,recv%domain,comm_icosa, mpi_req(ireq),ierr) |
---|
800 | ENDDO |
---|
801 | |
---|
802 | ENDDO |
---|
803 | |
---|
804 | CALL MPI_WAITALL(nreq,mpi_req,status,ierr) |
---|
805 | ! DO ind=1,ndomain |
---|
806 | ! field(ind)%rval2d(:)=0. |
---|
807 | ! ENDDO |
---|
808 | |
---|
809 | DO ind=1,ndomain |
---|
810 | rval4d=>field(ind)%rval4d |
---|
811 | |
---|
812 | req=>request(ind) |
---|
813 | DO isend=1,req%nsend |
---|
814 | send=>req%send(isend) |
---|
815 | DEALLOCATE(send%buffer_r4) |
---|
816 | ENDDO |
---|
817 | |
---|
818 | DO irecv=1,req%nrecv |
---|
819 | recv=>req%recv(irecv) |
---|
820 | buffer_r4=>recv%buffer_r4 |
---|
821 | value=>recv%value |
---|
822 | sgn=>recv%sign |
---|
823 | DO n=1,recv%size |
---|
824 | rval4d(value(n),:,:)=buffer_r4(n,:,:)*sgn(n) |
---|
825 | ENDDO |
---|
826 | DEALLOCATE(recv%buffer_r4) |
---|
827 | ENDDO |
---|
828 | |
---|
829 | ENDDO |
---|
830 | |
---|
831 | ENDIF |
---|
832 | |
---|
833 | ENDIF |
---|
834 | |
---|
835 | CALL trace_end("transfert_mpi") |
---|
836 | |
---|
837 | END SUBROUTINE transfert_request_mpi |
---|
838 | |
---|
839 | SUBROUTINE transfert_request(field,request) |
---|
840 | USE field_mod |
---|
841 | USE domain_mod |
---|
842 | IMPLICIT NONE |
---|
843 | TYPE(t_field),POINTER :: field(:) |
---|
844 | TYPE(t_request),POINTER :: request(:) |
---|
845 | REAL(rstd),POINTER :: rval2d(:) |
---|
846 | REAL(rstd),POINTER :: rval3d(:,:) |
---|
847 | REAL(rstd),POINTER :: rval4d(:,:,:) |
---|
848 | INTEGER :: ind |
---|
849 | TYPE(t_request),POINTER :: req |
---|
850 | INTEGER :: n |
---|
851 | REAL(rstd) :: var1,var2 |
---|
852 | |
---|
853 | DO ind=1,ndomain |
---|
854 | req=>request(ind) |
---|
855 | rval2d=>field(ind)%rval2d |
---|
856 | rval3d=>field(ind)%rval3d |
---|
857 | rval4d=>field(ind)%rval4d |
---|
858 | |
---|
859 | IF (field(ind)%data_type==type_real) THEN |
---|
860 | IF (field(ind)%ndim==2) THEN |
---|
861 | DO n=1,req%size |
---|
862 | rval2d(req%target_ind(n))=field(req%src_domain(n))%rval2d(req%src_ind(n))*req%target_sign(n) |
---|
863 | ENDDO |
---|
864 | ELSE IF (field(ind)%ndim==3) THEN |
---|
865 | DO n=1,req%size |
---|
866 | rval3d(req%target_ind(n),:)=field(req%src_domain(n))%rval3d(req%src_ind(n),:)*req%target_sign(n) |
---|
867 | ENDDO |
---|
868 | ELSE IF (field(ind)%ndim==4) THEN |
---|
869 | DO n=1,req%size |
---|
870 | rval4d(req%target_ind(n),:,:)=field(req%src_domain(n))%rval4d(req%src_ind(n),:,:)*req%target_sign(n) |
---|
871 | ENDDO |
---|
872 | ENDIF |
---|
873 | ENDIF |
---|
874 | |
---|
875 | ENDDO |
---|
876 | |
---|
877 | END SUBROUTINE transfert_request |
---|
878 | |
---|
879 | |
---|
880 | SUBROUTINE gather_field(field_loc,field_glo) |
---|
881 | USE field_mod |
---|
882 | USE domain_mod |
---|
883 | USE mpi_mod |
---|
884 | USE mpipara |
---|
885 | IMPLICIT NONE |
---|
886 | TYPE(t_field),POINTER :: field_loc(:) |
---|
887 | TYPE(t_field),POINTER :: field_glo(:) |
---|
888 | INTEGER, ALLOCATABLE :: mpi_req(:) |
---|
889 | INTEGER, ALLOCATABLE :: status(:,:) |
---|
890 | INTEGER :: ireq,nreq |
---|
891 | INTEGER :: ind_glo,ind_loc |
---|
892 | |
---|
893 | IF (.NOT. using_mpi) THEN |
---|
894 | |
---|
895 | DO ind_loc=1,ndomain |
---|
896 | IF (field_loc(ind_loc)%ndim==2) field_glo(ind_loc)%rval2d=field_loc(ind_loc)%rval2d |
---|
897 | IF (field_loc(ind_loc)%ndim==3) field_glo(ind_loc)%rval3d=field_loc(ind_loc)%rval3d |
---|
898 | IF (field_loc(ind_loc)%ndim==4) field_glo(ind_loc)%rval4d=field_loc(ind_loc)%rval4d |
---|
899 | ENDDO |
---|
900 | |
---|
901 | ELSE |
---|
902 | |
---|
903 | nreq=ndomain |
---|
904 | IF (mpi_rank==0) nreq=nreq+ndomain_glo |
---|
905 | ALLOCATE(mpi_req(nreq)) |
---|
906 | ALLOCATE(status(MPI_STATUS_SIZE,nreq)) |
---|
907 | |
---|
908 | |
---|
909 | ireq=0 |
---|
910 | IF (mpi_rank==0) THEN |
---|
911 | DO ind_glo=1,ndomain_glo |
---|
912 | ireq=ireq+1 |
---|
913 | |
---|
914 | IF (field_glo(ind_glo)%ndim==2) THEN |
---|
915 | CALL MPI_IRECV(field_glo(ind_glo)%rval2d,size(field_glo(ind_glo)%rval2d) , MPI_REAL8 , & |
---|
916 | domglo_rank(ind_glo),domglo_loc_ind(ind_glo), comm_icosa, mpi_req(ireq), ierr) |
---|
917 | |
---|
918 | ELSE IF (field_glo(ind_glo)%ndim==3) THEN |
---|
919 | CALL MPI_IRECV(field_glo(ind_glo)%rval3d,size(field_glo(ind_glo)%rval3d) , MPI_REAL8 , & |
---|
920 | domglo_rank(ind_glo),domglo_loc_ind(ind_glo), comm_icosa, mpi_req(ireq), ierr) |
---|
921 | |
---|
922 | ELSE IF (field_glo(ind_glo)%ndim==4) THEN |
---|
923 | CALL MPI_IRECV(field_glo(ind_glo)%rval4d,size(field_glo(ind_glo)%rval4d) , MPI_REAL8 , & |
---|
924 | domglo_rank(ind_glo),domglo_loc_ind(ind_glo), comm_icosa, mpi_req(ireq), ierr) |
---|
925 | ENDIF |
---|
926 | |
---|
927 | ENDDO |
---|
928 | ENDIF |
---|
929 | |
---|
930 | DO ind_loc=1,ndomain |
---|
931 | ireq=ireq+1 |
---|
932 | |
---|
933 | IF (field_loc(ind_loc)%ndim==2) THEN |
---|
934 | CALL MPI_ISEND(field_loc(ind_loc)%rval2d,size(field_loc(ind_loc)%rval2d) , MPI_REAL8 , & |
---|
935 | 0, ind_loc, comm_icosa, mpi_req(ireq), ierr) |
---|
936 | ELSE IF (field_loc(ind_loc)%ndim==3) THEN |
---|
937 | CALL MPI_ISEND(field_loc(ind_loc)%rval3d,size(field_loc(ind_loc)%rval3d) , MPI_REAL8 , & |
---|
938 | 0, ind_loc, comm_icosa, mpi_req(ireq), ierr) |
---|
939 | ELSE IF (field_loc(ind_loc)%ndim==4) THEN |
---|
940 | CALL MPI_ISEND(field_loc(ind_loc)%rval4d,size(field_loc(ind_loc)%rval4d) , MPI_REAL8 , & |
---|
941 | 0, ind_loc, comm_icosa, mpi_req(ireq), ierr) |
---|
942 | ENDIF |
---|
943 | |
---|
944 | ENDDO |
---|
945 | |
---|
946 | CALL MPI_WAITALL(nreq,mpi_req,status,ierr) |
---|
947 | |
---|
948 | ENDIF |
---|
949 | |
---|
950 | END SUBROUTINE gather_field |
---|
951 | |
---|
952 | |
---|
953 | END MODULE transfert_mpi_mod |
---|
954 | |
---|
955 | |
---|
956 | |
---|
957 | |
---|
958 | |
---|