1 | ! |
---|
2 | ! $Id$ |
---|
3 | ! |
---|
4 | ! AGRIF (Adaptive Grid Refinement In Fortran) |
---|
5 | ! |
---|
6 | ! Copyright (C) 2003 Laurent Debreu (Laurent.Debreu@imag.fr) |
---|
7 | ! Christophe Vouland (Christophe.Vouland@imag.fr) |
---|
8 | ! |
---|
9 | ! This program is free software; you can redistribute it and/or modify |
---|
10 | ! it under the terms of the GNU General Public License as published by |
---|
11 | ! the Free Software Foundation; either version 2 of the License, or |
---|
12 | ! (at your option) any later version. |
---|
13 | ! |
---|
14 | ! This program is distributed in the hope that it will be useful, |
---|
15 | ! but WITHOUT ANY WARRANTY; without even the implied warranty of |
---|
16 | ! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the |
---|
17 | ! GNU General Public License for more details. |
---|
18 | ! |
---|
19 | ! You should have received a copy of the GNU General Public License |
---|
20 | ! along with this program; if not, write to the Free Software |
---|
21 | ! Foundation, Inc., 59 Temple Place- Suite 330, Boston, MA 02111-1307, USA. |
---|
22 | ! |
---|
23 | ! |
---|
24 | module Agrif_Mpp |
---|
25 | ! |
---|
26 | use Agrif_Arrays |
---|
27 | use Agrif_Grids |
---|
28 | ! |
---|
29 | implicit none |
---|
30 | ! |
---|
31 | interface |
---|
32 | subroutine Agrif_get_proc_info ( imin, imax, jmin, jmax ) |
---|
33 | integer, intent(out) :: imin, imax |
---|
34 | integer, intent(out) :: jmin, jmax |
---|
35 | end subroutine Agrif_get_proc_info |
---|
36 | end interface |
---|
37 | ! |
---|
38 | integer, private :: Agrif_MPI_prec |
---|
39 | ! |
---|
40 | private :: Agrif_get_proc_info |
---|
41 | ! |
---|
42 | contains |
---|
43 | ! |
---|
44 | #if defined AGRIF_MPI |
---|
45 | !=================================================================================================== |
---|
46 | ! subroutine Agrif_MPI_Init |
---|
47 | !--------------------------------------------------------------------------------------------------- |
---|
48 | subroutine Agrif_MPI_Init ( comm ) |
---|
49 | !--------------------------------------------------------------------------------------------------- |
---|
50 | integer, optional, intent(in) :: comm !< MPI communicator to be attached to the root grid. |
---|
51 | ! |
---|
52 | include 'mpif.h' |
---|
53 | integer :: code, ierr |
---|
54 | logical :: mpi_was_called |
---|
55 | integer :: current_mpi_prec |
---|
56 | ! |
---|
57 | call MPI_INITIALIZED( mpi_was_called, code ) |
---|
58 | if( code /= MPI_SUCCESS ) then |
---|
59 | write(*,*) ': Error in routine mpi_initialized' |
---|
60 | call MPI_ABORT( MPI_COMM_WORLD, code, ierr ) |
---|
61 | endif |
---|
62 | if( .not. mpi_was_called ) then |
---|
63 | write(*,*) '### AGRIF Error : you should call Agrif_MPI_Init *after* MPI_Init.' |
---|
64 | stop |
---|
65 | endif |
---|
66 | |
---|
67 | current_mpi_prec = KIND(1.0) |
---|
68 | if (current_mpi_prec == 4) then |
---|
69 | Agrif_MPI_prec = MPI_REAL4 |
---|
70 | else |
---|
71 | Agrif_MPI_prec = MPI_REAL8 |
---|
72 | endif |
---|
73 | ! |
---|
74 | if ( present(comm) ) then |
---|
75 | call Agrif_MPI_switch_comm(comm) |
---|
76 | else |
---|
77 | call Agrif_MPI_switch_comm(MPI_COMM_WORLD) |
---|
78 | endif |
---|
79 | ! |
---|
80 | Agrif_Mygrid % communicator = Agrif_mpi_comm |
---|
81 | ! |
---|
82 | if ( Agrif_Parallel_sisters ) then |
---|
83 | call Agrif_Init_ProcList( Agrif_Mygrid % proc_def_list, Agrif_Nbprocs ) |
---|
84 | call Agrif_pl_copy( Agrif_Mygrid % proc_def_list, Agrif_Mygrid % required_proc_list ) |
---|
85 | endif |
---|
86 | !--------------------------------------------------------------------------------------------------- |
---|
87 | end subroutine Agrif_MPI_Init |
---|
88 | !=================================================================================================== |
---|
89 | ! |
---|
90 | !=================================================================================================== |
---|
91 | subroutine Agrif_MPI_switch_comm ( comm ) |
---|
92 | !--------------------------------------------------------------------------------------------------- |
---|
93 | integer, intent(in) :: comm !< MPI communicator you want to switch to. |
---|
94 | ! |
---|
95 | include 'mpif.h' |
---|
96 | integer :: code |
---|
97 | logical :: mpi_was_called |
---|
98 | ! |
---|
99 | call MPI_INITIALIZED( mpi_was_called, code ) |
---|
100 | if ( .not. mpi_was_called ) return |
---|
101 | ! |
---|
102 | call MPI_COMM_SIZE(comm, Agrif_Nbprocs, code) |
---|
103 | call MPI_COMM_RANK(comm, Agrif_ProcRank, code) |
---|
104 | Agrif_mpi_comm = comm |
---|
105 | !--------------------------------------------------------------------------------------------------- |
---|
106 | end subroutine Agrif_MPI_switch_comm |
---|
107 | !=================================================================================================== |
---|
108 | ! |
---|
109 | !=================================================================================================== |
---|
110 | function Agrif_MPI_get_grid_comm ( ) result ( comm ) |
---|
111 | !--------------------------------------------------------------------------------------------------- |
---|
112 | integer :: comm |
---|
113 | comm = Agrif_Curgrid % communicator |
---|
114 | !--------------------------------------------------------------------------------------------------- |
---|
115 | end function Agrif_MPI_get_grid_comm |
---|
116 | !=================================================================================================== |
---|
117 | ! |
---|
118 | !=================================================================================================== |
---|
119 | subroutine Agrif_MPI_set_grid_comm ( comm ) |
---|
120 | !--------------------------------------------------------------------------------------------------- |
---|
121 | integer, intent(in) :: comm |
---|
122 | Agrif_Curgrid % communicator = comm |
---|
123 | !--------------------------------------------------------------------------------------------------- |
---|
124 | end subroutine Agrif_MPI_set_grid_comm |
---|
125 | !=================================================================================================== |
---|
126 | ! |
---|
127 | !=================================================================================================== |
---|
128 | subroutine Agrif_Init_ProcList ( proclist, nbprocs ) |
---|
129 | !--------------------------------------------------------------------------------------------------- |
---|
130 | type(Agrif_Proc_List), intent(inout) :: proclist |
---|
131 | integer, intent(in) :: nbprocs |
---|
132 | ! |
---|
133 | include 'mpif.h' |
---|
134 | type(Agrif_Proc), pointer :: new_proc |
---|
135 | integer :: p, ierr |
---|
136 | integer :: imin, imax, jmin, jmax |
---|
137 | integer, dimension(5) :: local_proc_grid_info |
---|
138 | integer, dimension(5,nbprocs) :: all_procs_grid_info |
---|
139 | ! |
---|
140 | call Agrif_get_proc_info(imin, imax, jmin, jmax) |
---|
141 | ! |
---|
142 | local_proc_grid_info(:) = (/Agrif_Procrank, imin, jmin, imax, jmax/) |
---|
143 | ! |
---|
144 | call MPI_ALLGATHER(local_proc_grid_info, 5, MPI_INTEGER, & |
---|
145 | all_procs_grid_info, 5, MPI_INTEGER, Agrif_mpi_comm, ierr) |
---|
146 | ! |
---|
147 | do p = 1,nbprocs |
---|
148 | ! |
---|
149 | allocate(new_proc) |
---|
150 | new_proc % pn = all_procs_grid_info(1,p) |
---|
151 | new_proc % imin(1) = all_procs_grid_info(2,p) |
---|
152 | new_proc % imin(2) = all_procs_grid_info(3,p) |
---|
153 | new_proc % imax(1) = all_procs_grid_info(4,p) |
---|
154 | new_proc % imax(2) = all_procs_grid_info(5,p) |
---|
155 | call Agrif_pl_append( proclist, new_proc ) |
---|
156 | ! |
---|
157 | enddo |
---|
158 | ! |
---|
159 | !--------------------------------------------------------------------------------------------------- |
---|
160 | end subroutine Agrif_Init_ProcList |
---|
161 | !=================================================================================================== |
---|
162 | ! |
---|
163 | !=================================================================================================== |
---|
164 | ! subroutine Get_External_Data_first |
---|
165 | !--------------------------------------------------------------------------------------------------- |
---|
166 | subroutine Get_External_Data_first ( pttruetab, cetruetab, pttruetabwhole, cetruetabwhole, & |
---|
167 | nbdim, memberoutall, coords, sendtoproc, recvfromproc, & |
---|
168 | imin, imax, imin_recv, imax_recv ) |
---|
169 | !--------------------------------------------------------------------------------------------------- |
---|
170 | include 'mpif.h' |
---|
171 | ! |
---|
172 | integer, intent(in) :: nbdim |
---|
173 | integer, dimension(nbdim,0:Agrif_NbProcs-1), intent(in) :: pttruetab, cetruetab |
---|
174 | integer, dimension(nbdim,0:Agrif_NbProcs-1), intent(in) :: pttruetabwhole,cetruetabwhole |
---|
175 | logical, dimension(0:Agrif_Nbprocs-1), intent(in) :: memberoutall |
---|
176 | integer, dimension(nbdim), intent(in) :: coords |
---|
177 | logical, dimension(0:Agrif_Nbprocs-1), intent(out) :: sendtoproc |
---|
178 | logical, dimension(0:Agrif_Nbprocs-1), intent(out) :: recvfromproc |
---|
179 | integer, dimension(nbdim,0:Agrif_NbProcs-1), intent(out) :: imin,imax |
---|
180 | integer, dimension(nbdim,0:Agrif_NbProcs-1), intent(out) :: imin_recv,imax_recv |
---|
181 | ! |
---|
182 | integer :: imintmp, imaxtmp, i, j, k, i1 |
---|
183 | integer :: imin1,imax1 |
---|
184 | logical :: tochange,tochangebis |
---|
185 | integer, dimension(nbdim,0:Agrif_NbProcs-1) :: pttruetab2,cetruetab2 |
---|
186 | ! |
---|
187 | ! pttruetab2 and cetruetab2 are modified arrays in order to always |
---|
188 | ! send the most inner points |
---|
189 | ! |
---|
190 | pttruetab2(:,Agrif_Procrank) = pttruetab(:,Agrif_Procrank) |
---|
191 | cetruetab2(:,Agrif_Procrank) = cetruetab(:,Agrif_Procrank) |
---|
192 | ! |
---|
193 | do k = 0,Agrif_Nbprocs-1 |
---|
194 | do i = 1,nbdim |
---|
195 | tochangebis = .TRUE. |
---|
196 | DO i1 = 1,nbdim |
---|
197 | IF (i /= i1) THEN |
---|
198 | IF ( (pttruetab(i1,Agrif_Procrank) /= pttruetab(i1,k)) .OR. & |
---|
199 | (cetruetab(i1,Agrif_Procrank) /= cetruetab(i1,k))) THEN |
---|
200 | tochangebis = .FALSE. |
---|
201 | EXIT |
---|
202 | ENDIF |
---|
203 | ENDIF |
---|
204 | ENDDO |
---|
205 | IF (tochangebis) THEN |
---|
206 | imin1 = max(pttruetab(i,Agrif_Procrank), pttruetab(i,k)) |
---|
207 | imax1 = min(cetruetab(i,Agrif_Procrank), cetruetab(i,k)) |
---|
208 | ! Always send the most interior points |
---|
209 | |
---|
210 | tochange = .false. |
---|
211 | IF (cetruetab(i,Agrif_Procrank) > cetruetab(i,k)) THEN |
---|
212 | DO j=imin1,imax1 |
---|
213 | IF ((cetruetab(i,k)-j) > (j-pttruetab(i,Agrif_Procrank))) THEN |
---|
214 | imintmp = j+1 |
---|
215 | tochange = .TRUE. |
---|
216 | ELSE |
---|
217 | EXIT |
---|
218 | ENDIF |
---|
219 | ENDDO |
---|
220 | ENDIF |
---|
221 | |
---|
222 | if (tochange) then |
---|
223 | pttruetab2(i,Agrif_Procrank) = imintmp |
---|
224 | endif |
---|
225 | |
---|
226 | tochange = .FALSE. |
---|
227 | imaxtmp=0 |
---|
228 | IF (pttruetab(i,Agrif_Procrank) < pttruetab(i,k)) THEN |
---|
229 | DO j=imax1,imin1,-1 |
---|
230 | IF ((j-pttruetab(i,k)) > (cetruetab(i,Agrif_Procrank)-j)) THEN |
---|
231 | imaxtmp = j-1 |
---|
232 | tochange = .TRUE. |
---|
233 | ELSE |
---|
234 | EXIT |
---|
235 | ENDIF |
---|
236 | ENDDO |
---|
237 | ENDIF |
---|
238 | |
---|
239 | if (tochange) then |
---|
240 | cetruetab2(i,Agrif_Procrank) = imaxtmp |
---|
241 | endif |
---|
242 | ENDIF |
---|
243 | enddo |
---|
244 | enddo |
---|
245 | |
---|
246 | do k = 0,Agrif_NbProcs-1 |
---|
247 | ! |
---|
248 | sendtoproc(k) = .true. |
---|
249 | ! |
---|
250 | !CDIR SHORTLOOP |
---|
251 | do i = 1,nbdim |
---|
252 | imin(i,k) = max(pttruetab2(i,Agrif_Procrank), pttruetabwhole(i,k)) |
---|
253 | imax(i,k) = min(cetruetab2(i,Agrif_Procrank), cetruetabwhole(i,k)) |
---|
254 | ! |
---|
255 | if ( (imin(i,k) > imax(i,k)) .and. (coords(i) /= 0) ) then |
---|
256 | sendtoproc(k) = .false. |
---|
257 | endif |
---|
258 | enddo |
---|
259 | IF ( .not. memberoutall(k) ) THEN |
---|
260 | sendtoproc(k) = .false. |
---|
261 | ENDIF |
---|
262 | enddo |
---|
263 | ! |
---|
264 | call Exchangesamelevel_first(sendtoproc,nbdim,imin,imax,recvfromproc,imin_recv,imax_recv) |
---|
265 | !--------------------------------------------------------------------------------------------------- |
---|
266 | end subroutine Get_External_Data_first |
---|
267 | !=================================================================================================== |
---|
268 | ! |
---|
269 | !=================================================================================================== |
---|
270 | ! subroutine ExchangeSameLevel_first |
---|
271 | !--------------------------------------------------------------------------------------------------- |
---|
272 | subroutine ExchangeSameLevel_first ( sendtoproc, nbdim, imin, imax, recvfromproc, & |
---|
273 | imin_recv, imax_recv ) |
---|
274 | !--------------------------------------------------------------------------------------------------- |
---|
275 | LOGICAL, DIMENSION(0:Agrif_Nbprocs-1), intent(in) :: sendtoproc |
---|
276 | INTEGER, intent(in) :: nbdim |
---|
277 | INTEGER, DIMENSION(nbdim,0:Agrif_Nbprocs-1), intent(in) :: imin |
---|
278 | INTEGER, DIMENSION(nbdim,0:Agrif_Nbprocs-1), intent(in) :: imax |
---|
279 | LOGICAL, DIMENSION(0:Agrif_Nbprocs-1), intent(out) :: recvfromproc |
---|
280 | INTEGER, DIMENSION(nbdim,0:Agrif_Nbprocs-1), intent(out) :: imin_recv |
---|
281 | INTEGER, DIMENSION(nbdim,0:Agrif_Nbprocs-1), intent(out) :: imax_recv |
---|
282 | ! |
---|
283 | include 'mpif.h' |
---|
284 | INTEGER :: k |
---|
285 | INTEGER :: etiquette = 100 |
---|
286 | INTEGER :: code |
---|
287 | LOGICAL :: res |
---|
288 | INTEGER, DIMENSION(MPI_STATUS_SIZE) :: statut |
---|
289 | INTEGER, DIMENSION(nbdim,2,0:Agrif_Nbprocs-1) :: iminmax_temp |
---|
290 | |
---|
291 | do k = 0,Agrif_ProcRank-1 |
---|
292 | ! |
---|
293 | call MPI_SEND(sendtoproc(k),1,MPI_LOGICAL,k,etiquette,Agrif_mpi_comm,code) |
---|
294 | ! |
---|
295 | if (sendtoproc(k)) then |
---|
296 | iminmax_temp(:,1,k) = imin(:,k) |
---|
297 | iminmax_temp(:,2,k) = imax(:,k) |
---|
298 | call MPI_SEND(iminmax_temp(:,:,k),2*nbdim,MPI_INTEGER,k,etiquette,Agrif_mpi_comm,code) |
---|
299 | endif |
---|
300 | ! |
---|
301 | enddo |
---|
302 | ! |
---|
303 | ! Reception from others processors of the necessary part of the parent grid |
---|
304 | do k = Agrif_ProcRank+1,Agrif_Nbprocs-1 |
---|
305 | ! |
---|
306 | call MPI_RECV(res,1,MPI_LOGICAL,k,etiquette,Agrif_mpi_comm,statut,code) |
---|
307 | recvfromproc(k) = res |
---|
308 | ! |
---|
309 | if (recvfromproc(k)) then |
---|
310 | call MPI_RECV(iminmax_temp(:,:,k),2*nbdim,MPI_INTEGER,k,etiquette, & |
---|
311 | Agrif_mpi_comm,statut,code) |
---|
312 | imin_recv(:,k) = iminmax_temp(:,1,k) |
---|
313 | imax_recv(:,k) = iminmax_temp(:,2,k) |
---|
314 | endif |
---|
315 | ! |
---|
316 | enddo |
---|
317 | |
---|
318 | ! Reception from others processors of the necessary part of the parent grid |
---|
319 | do k = Agrif_ProcRank+1,Agrif_Nbprocs-1 |
---|
320 | ! |
---|
321 | call MPI_SEND(sendtoproc(k),1,MPI_LOGICAL,k,etiquette,Agrif_mpi_comm,code) |
---|
322 | ! |
---|
323 | if (sendtoproc(k)) then |
---|
324 | ! |
---|
325 | iminmax_temp(:,1,k) = imin(:,k) |
---|
326 | iminmax_temp(:,2,k) = imax(:,k) |
---|
327 | |
---|
328 | call MPI_SEND(iminmax_temp(:,:,k),2*nbdim,MPI_INTEGER,k,etiquette, & |
---|
329 | Agrif_mpi_comm,code) |
---|
330 | endif |
---|
331 | ! |
---|
332 | enddo |
---|
333 | ! |
---|
334 | ! |
---|
335 | ! Reception from others processors of the necessary part of the parent grid |
---|
336 | do k = Agrif_ProcRank-1,0,-1 |
---|
337 | ! |
---|
338 | call MPI_RECV(res,1,MPI_LOGICAL,k,etiquette,Agrif_mpi_comm,statut,code) |
---|
339 | recvfromproc(k) = res |
---|
340 | ! |
---|
341 | if (recvfromproc(k)) then |
---|
342 | ! |
---|
343 | call MPI_RECV(iminmax_temp(:,:,k),2*nbdim,MPI_INTEGER,k,etiquette, & |
---|
344 | Agrif_mpi_comm,statut,code) |
---|
345 | |
---|
346 | imin_recv(:,k) = iminmax_temp(:,1,k) |
---|
347 | imax_recv(:,k) = iminmax_temp(:,2,k) |
---|
348 | endif |
---|
349 | ! |
---|
350 | enddo |
---|
351 | !--------------------------------------------------------------------------------------------------- |
---|
352 | end subroutine ExchangeSamelevel_first |
---|
353 | !=================================================================================================== |
---|
354 | ! |
---|
355 | !=================================================================================================== |
---|
356 | ! subroutine ExchangeSameLevel |
---|
357 | !--------------------------------------------------------------------------------------------------- |
---|
358 | subroutine ExchangeSameLevel ( sendtoproc, recvfromproc, nbdim, & |
---|
359 | pttruetabwhole, cetruetabwhole, & |
---|
360 | imin, imax, imin_recv, imax_recv, & |
---|
361 | memberout, tempC, tempCextend ) |
---|
362 | !--------------------------------------------------------------------------------------------------- |
---|
363 | LOGICAL, DIMENSION(0:Agrif_Nbprocs-1), intent(in) :: sendtoproc |
---|
364 | LOGICAL, DIMENSION(0:Agrif_Nbprocs-1), intent(in) :: recvfromproc |
---|
365 | INTEGER, intent(in) :: nbdim |
---|
366 | INTEGER, DIMENSION(nbdim,0:Agrif_Nbprocs-1), intent(in) :: pttruetabwhole |
---|
367 | INTEGER, DIMENSION(nbdim,0:Agrif_Nbprocs-1), intent(in) :: cetruetabwhole |
---|
368 | INTEGER, DIMENSION(nbdim,0:Agrif_Nbprocs-1), intent(in) :: imin, imax |
---|
369 | INTEGER, DIMENSION(nbdim,0:Agrif_Nbprocs-1), intent(in) :: imin_recv, imax_recv |
---|
370 | LOGICAL, intent(in) :: memberout |
---|
371 | TYPE(Agrif_Variable), pointer, intent(inout) :: tempC, tempCextend |
---|
372 | ! |
---|
373 | include 'mpif.h' |
---|
374 | INTEGER :: i,k |
---|
375 | INTEGER :: etiquette = 100 |
---|
376 | INTEGER :: code, datasize |
---|
377 | INTEGER, DIMENSION(MPI_STATUS_SIZE) :: statut |
---|
378 | TYPE(Agrif_Variable), pointer, SAVE :: temprecv |
---|
379 | ! |
---|
380 | IF (memberout) THEN |
---|
381 | call Agrif_array_allocate(tempCextend, pttruetabwhole(:,Agrif_ProcRank), & |
---|
382 | cetruetabwhole(:,Agrif_ProcRank),nbdim) |
---|
383 | call Agrif_var_set_array_tozero(tempCextend,nbdim) |
---|
384 | ENDIF |
---|
385 | ! |
---|
386 | IF (sendtoproc(Agrif_ProcRank)) THEN |
---|
387 | call Agrif_var_copy_array(tempCextend,imin(:,Agrif_Procrank),imax(:,Agrif_Procrank), & |
---|
388 | tempC, imin(:,Agrif_Procrank),imax(:,Agrif_Procrank), & |
---|
389 | nbdim) |
---|
390 | ENDIF |
---|
391 | ! |
---|
392 | do k = 0,Agrif_ProcRank-1 |
---|
393 | ! |
---|
394 | if (sendtoproc(k)) then |
---|
395 | ! |
---|
396 | datasize = 1 |
---|
397 | ! |
---|
398 | !CDIR SHORTLOOP |
---|
399 | do i = 1,nbdim |
---|
400 | datasize = datasize * (imax(i,k)-imin(i,k)+1) |
---|
401 | enddo |
---|
402 | ! |
---|
403 | SELECT CASE(nbdim) |
---|
404 | CASE(1) |
---|
405 | call MPI_SEND(tempC%array1(imin(1,k):imax(1,k)), & |
---|
406 | datasize,Agrif_MPI_prec,k,etiquette, & |
---|
407 | Agrif_mpi_comm,code) |
---|
408 | CASE(2) |
---|
409 | call MPI_SEND(tempC%array2(imin(1,k):imax(1,k), & |
---|
410 | imin(2,k):imax(2,k)), & |
---|
411 | datasize,Agrif_MPI_prec,k,etiquette, & |
---|
412 | Agrif_mpi_comm,code) |
---|
413 | CASE(3) |
---|
414 | call Agrif_Send_3Darray(tempC%array3,lbound(tempC%array3),imin(:,k),imax(:,k),k) |
---|
415 | CASE(4) |
---|
416 | call MPI_SEND(tempC%array4(imin(1,k):imax(1,k), & |
---|
417 | imin(2,k):imax(2,k), & |
---|
418 | imin(3,k):imax(3,k), & |
---|
419 | imin(4,k):imax(4,k)), & |
---|
420 | datasize,Agrif_MPI_prec,k,etiquette, & |
---|
421 | Agrif_mpi_comm,code) |
---|
422 | CASE(5) |
---|
423 | call MPI_SEND(tempC%array5(imin(1,k):imax(1,k), & |
---|
424 | imin(2,k):imax(2,k), & |
---|
425 | imin(3,k):imax(3,k), & |
---|
426 | imin(4,k):imax(4,k), & |
---|
427 | imin(5,k):imax(5,k)), & |
---|
428 | datasize,Agrif_MPI_prec,k,etiquette, & |
---|
429 | Agrif_mpi_comm,code) |
---|
430 | CASE(6) |
---|
431 | call MPI_SEND(tempC%array6(imin(1,k):imax(1,k), & |
---|
432 | imin(2,k):imax(2,k), & |
---|
433 | imin(3,k):imax(3,k), & |
---|
434 | imin(4,k):imax(4,k), & |
---|
435 | imin(5,k):imax(5,k), & |
---|
436 | imin(6,k):imax(6,k)), & |
---|
437 | datasize,Agrif_MPI_prec,k,etiquette, & |
---|
438 | Agrif_mpi_comm,code) |
---|
439 | END SELECT |
---|
440 | ! |
---|
441 | endif |
---|
442 | enddo |
---|
443 | ! |
---|
444 | ! Reception from others processors of the necessary part of the parent grid |
---|
445 | do k = Agrif_ProcRank+1,Agrif_Nbprocs-1 |
---|
446 | ! |
---|
447 | if (recvfromproc(k)) then |
---|
448 | ! |
---|
449 | datasize = 1 |
---|
450 | ! |
---|
451 | !CDIR SHORTLOOP |
---|
452 | do i = 1,nbdim |
---|
453 | datasize = datasize * (imax_recv(i,k)-imin_recv(i,k)+1) |
---|
454 | enddo |
---|
455 | |
---|
456 | if (.not.associated(temprecv)) allocate(temprecv) |
---|
457 | call Agrif_array_allocate(temprecv,imin_recv(:,k),imax_recv(:,k),nbdim) |
---|
458 | |
---|
459 | SELECT CASE(nbdim) |
---|
460 | CASE(1) |
---|
461 | call MPI_RECV(temprecv%array1,datasize,Agrif_MPI_prec,k,etiquette, & |
---|
462 | Agrif_mpi_comm,statut,code) |
---|
463 | CASE(2) |
---|
464 | call MPI_RECV(temprecv%array2,datasize,Agrif_MPI_prec,k,etiquette, & |
---|
465 | Agrif_mpi_comm,statut,code) |
---|
466 | CASE(3) |
---|
467 | call MPI_RECV(temprecv%array3,datasize,Agrif_MPI_prec,k,etiquette, & |
---|
468 | Agrif_mpi_comm,statut,code) |
---|
469 | CASE(4) |
---|
470 | call MPI_RECV(temprecv%array4,datasize,Agrif_MPI_prec,k,etiquette, & |
---|
471 | Agrif_mpi_comm,statut,code) |
---|
472 | CASE(5) |
---|
473 | call MPI_RECV(temprecv%array5,datasize,Agrif_MPI_prec,k,etiquette, & |
---|
474 | Agrif_mpi_comm,statut,code) |
---|
475 | CASE(6) |
---|
476 | call MPI_RECV(temprecv%array6,datasize,Agrif_MPI_prec,k,etiquette, & |
---|
477 | Agrif_mpi_comm,statut,code) |
---|
478 | END SELECT |
---|
479 | |
---|
480 | call Agrif_var_replace_value(tempCextend,temprecv,imin_recv(:,k),imax_recv(:,k),0.,nbdim) |
---|
481 | call Agrif_array_deallocate(temprecv,nbdim) |
---|
482 | ! |
---|
483 | endif |
---|
484 | enddo |
---|
485 | |
---|
486 | ! Reception from others processors of the necessary part of the parent grid |
---|
487 | do k = Agrif_ProcRank+1,Agrif_Nbprocs-1 |
---|
488 | ! |
---|
489 | if (sendtoproc(k)) then |
---|
490 | ! |
---|
491 | SELECT CASE(nbdim) |
---|
492 | CASE(1) |
---|
493 | datasize=SIZE(tempC%array1(imin(1,k):imax(1,k))) |
---|
494 | call MPI_SEND(tempC%array1(imin(1,k):imax(1,k)), & |
---|
495 | datasize,Agrif_MPI_prec,k,etiquette, & |
---|
496 | Agrif_mpi_comm,code) |
---|
497 | CASE(2) |
---|
498 | datasize=SIZE(tempC%array2(imin(1,k):imax(1,k), & |
---|
499 | imin(2,k):imax(2,k))) |
---|
500 | call MPI_SEND(tempC%array2(imin(1,k):imax(1,k), & |
---|
501 | imin(2,k):imax(2,k)), & |
---|
502 | datasize,Agrif_MPI_prec,k,etiquette, & |
---|
503 | Agrif_mpi_comm,code) |
---|
504 | CASE(3) |
---|
505 | datasize=SIZE(tempC%array3(imin(1,k):imax(1,k), & |
---|
506 | imin(2,k):imax(2,k), & |
---|
507 | imin(3,k):imax(3,k))) |
---|
508 | call MPI_SEND(tempC%array3(imin(1,k):imax(1,k), & |
---|
509 | imin(2,k):imax(2,k), & |
---|
510 | imin(3,k):imax(3,k)), & |
---|
511 | datasize,Agrif_MPI_prec,k,etiquette, & |
---|
512 | Agrif_mpi_comm,code) |
---|
513 | CASE(4) |
---|
514 | datasize=SIZE(tempC%array4(imin(1,k):imax(1,k), & |
---|
515 | imin(2,k):imax(2,k), & |
---|
516 | imin(3,k):imax(3,k), & |
---|
517 | imin(4,k):imax(4,k))) |
---|
518 | call MPI_SEND(tempC%array4(imin(1,k):imax(1,k), & |
---|
519 | imin(2,k):imax(2,k), & |
---|
520 | imin(3,k):imax(3,k), & |
---|
521 | imin(4,k):imax(4,k)), & |
---|
522 | datasize,Agrif_MPI_prec,k,etiquette, & |
---|
523 | Agrif_mpi_comm,code) |
---|
524 | CASE(5) |
---|
525 | datasize=SIZE(tempC%array5(imin(1,k):imax(1,k), & |
---|
526 | imin(2,k):imax(2,k), & |
---|
527 | imin(3,k):imax(3,k), & |
---|
528 | imin(4,k):imax(4,k), & |
---|
529 | imin(5,k):imax(5,k))) |
---|
530 | call MPI_SEND(tempC%array5(imin(1,k):imax(1,k), & |
---|
531 | imin(2,k):imax(2,k), & |
---|
532 | imin(3,k):imax(3,k), & |
---|
533 | imin(4,k):imax(4,k), & |
---|
534 | imin(5,k):imax(5,k)), & |
---|
535 | datasize,Agrif_MPI_prec,k,etiquette, & |
---|
536 | Agrif_mpi_comm,code) |
---|
537 | CASE(6) |
---|
538 | datasize=SIZE(tempC%array6(imin(1,k):imax(1,k), & |
---|
539 | imin(2,k):imax(2,k), & |
---|
540 | imin(3,k):imax(3,k), & |
---|
541 | imin(4,k):imax(4,k), & |
---|
542 | imin(5,k):imax(5,k), & |
---|
543 | imin(6,k):imax(6,k))) |
---|
544 | call MPI_SEND(tempC%array6(imin(1,k):imax(1,k), & |
---|
545 | imin(2,k):imax(2,k), & |
---|
546 | imin(3,k):imax(3,k), & |
---|
547 | imin(4,k):imax(4,k), & |
---|
548 | imin(5,k):imax(5,k), & |
---|
549 | imin(6,k):imax(6,k)), & |
---|
550 | datasize,Agrif_MPI_prec,k,etiquette, & |
---|
551 | Agrif_mpi_comm,code) |
---|
552 | END SELECT |
---|
553 | ! |
---|
554 | endif |
---|
555 | ! |
---|
556 | enddo |
---|
557 | ! |
---|
558 | ! Reception from others processors of the necessary part of the parent grid |
---|
559 | do k = Agrif_ProcRank-1,0,-1 |
---|
560 | ! |
---|
561 | if (recvfromproc(k)) then |
---|
562 | ! |
---|
563 | if (.not.associated(temprecv)) allocate(temprecv) |
---|
564 | call Agrif_array_allocate(temprecv,imin_recv(:,k),imax_recv(:,k),nbdim) |
---|
565 | |
---|
566 | SELECT CASE(nbdim) |
---|
567 | CASE(1) |
---|
568 | datasize=SIZE(temprecv%array1) |
---|
569 | call MPI_RECV(temprecv%array1,datasize,Agrif_MPI_prec,k,etiquette,& |
---|
570 | Agrif_mpi_comm,statut,code) |
---|
571 | CASE(2) |
---|
572 | datasize=SIZE(temprecv%array2) |
---|
573 | call MPI_RECV(temprecv%array2,datasize,Agrif_MPI_prec,k,etiquette,& |
---|
574 | Agrif_mpi_comm,statut,code) |
---|
575 | CASE(3) |
---|
576 | datasize=SIZE(temprecv%array3) |
---|
577 | call MPI_RECV(temprecv%array3,datasize,Agrif_MPI_prec,k,etiquette,& |
---|
578 | Agrif_mpi_comm,statut,code) |
---|
579 | CASE(4) |
---|
580 | datasize=SIZE(temprecv%array4) |
---|
581 | call MPI_RECV(temprecv%array4,datasize,Agrif_MPI_prec,k,etiquette,& |
---|
582 | Agrif_mpi_comm,statut,code) |
---|
583 | CASE(5) |
---|
584 | datasize=SIZE(temprecv%array5) |
---|
585 | call MPI_RECV(temprecv%array5,datasize,Agrif_MPI_prec,k,etiquette,& |
---|
586 | Agrif_mpi_comm,statut,code) |
---|
587 | CASE(6) |
---|
588 | datasize=SIZE(temprecv%array6) |
---|
589 | call MPI_RECV(temprecv%array6,datasize,Agrif_MPI_prec,k,etiquette,& |
---|
590 | Agrif_mpi_comm,statut,code) |
---|
591 | END SELECT |
---|
592 | |
---|
593 | call Agrif_var_replace_value(tempCextend,temprecv,imin_recv(:,k),imax_recv(:,k),0.,nbdim) |
---|
594 | call Agrif_array_deallocate(temprecv,nbdim) |
---|
595 | ! |
---|
596 | endif |
---|
597 | ! |
---|
598 | enddo |
---|
599 | !--------------------------------------------------------------------------------------------------- |
---|
600 | end subroutine ExchangeSamelevel |
---|
601 | !=================================================================================================== |
---|
602 | ! |
---|
603 | !=================================================================================================== |
---|
604 | ! subroutine Agrif_Send_3Darray |
---|
605 | !--------------------------------------------------------------------------------------------------- |
---|
606 | subroutine Agrif_Send_3Darray ( tab3D, bounds, imin, imax, k ) |
---|
607 | !--------------------------------------------------------------------------------------------------- |
---|
608 | integer, dimension(3), intent(in) :: bounds |
---|
609 | real, dimension(bounds(1):,bounds(2):,bounds(3):), target, intent(in) :: tab3D |
---|
610 | integer, dimension(3), intent(in) :: imin, imax |
---|
611 | integer, intent(in) :: k |
---|
612 | ! |
---|
613 | integer :: etiquette = 100 |
---|
614 | integer :: datasize, code |
---|
615 | include 'mpif.h' |
---|
616 | |
---|
617 | datasize = SIZE(tab3D(imin(1):imax(1), & |
---|
618 | imin(2):imax(2), & |
---|
619 | imin(3):imax(3))) |
---|
620 | |
---|
621 | call MPI_SEND( tab3D( imin(1):imax(1), & |
---|
622 | imin(2):imax(2), & |
---|
623 | imin(3):imax(3)), & |
---|
624 | datasize,Agrif_MPI_prec,k,etiquette,Agrif_mpi_comm,code) |
---|
625 | !--------------------------------------------------------------------------------------------------- |
---|
626 | end subroutine Agrif_Send_3Darray |
---|
627 | !=================================================================================================== |
---|
628 | ! |
---|
629 | #else |
---|
630 | subroutine dummy_Agrif_Mpp () |
---|
631 | end subroutine dummy_Agrif_Mpp |
---|
632 | #endif |
---|
633 | ! |
---|
634 | end Module Agrif_Mpp |
---|