[4777] | 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, & |
---|
[13027] | 168 | imin, imax, imin_recv, imax_recv, bornesmin, bornesmax ) |
---|
[4777] | 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 |
---|
[13027] | 181 | integer, dimension(nbdim,0:Agrif_NbProcs-1), intent(in) :: bornesmin, bornesmax |
---|
[4777] | 182 | ! |
---|
| 183 | integer :: imintmp, imaxtmp, i, j, k, i1 |
---|
| 184 | integer :: imin1,imax1 |
---|
| 185 | logical :: tochange,tochangebis |
---|
| 186 | integer, dimension(nbdim,0:Agrif_NbProcs-1) :: pttruetab2,cetruetab2 |
---|
| 187 | ! |
---|
| 188 | ! pttruetab2 and cetruetab2 are modified arrays in order to always |
---|
| 189 | ! send the most inner points |
---|
| 190 | ! |
---|
| 191 | pttruetab2(:,Agrif_Procrank) = pttruetab(:,Agrif_Procrank) |
---|
| 192 | cetruetab2(:,Agrif_Procrank) = cetruetab(:,Agrif_Procrank) |
---|
[13027] | 193 | |
---|
| 194 | if (agrif_debug_interp) then |
---|
| 195 | print *,'DANS Get_External_Data_first avec proc : ',Agrif_Procrank |
---|
| 196 | do k=0,Agrif_Nbprocs-1 |
---|
| 197 | print *,'Processeur ',k |
---|
| 198 | do i=1,nbdim |
---|
| 199 | print *,'ptcetretab = ',i,pttruetab(i,k),cetruetab(i,k) |
---|
| 200 | print *,'ptcetretabwhole = ',i,pttruetabwhole(i,k),cetruetabwhole(i,k) |
---|
| 201 | enddo |
---|
| 202 | enddo |
---|
| 203 | endif |
---|
[4777] | 204 | ! |
---|
| 205 | do k = 0,Agrif_Nbprocs-1 |
---|
[13027] | 206 | if (agrif_debug_interp) then |
---|
| 207 | print *,'Proc : ',k |
---|
| 208 | endif |
---|
[4777] | 209 | do i = 1,nbdim |
---|
[13027] | 210 | if (agrif_debug_interp) then |
---|
| 211 | print *,'Direction : ',i |
---|
| 212 | endif |
---|
[4777] | 213 | tochangebis = .TRUE. |
---|
| 214 | DO i1 = 1,nbdim |
---|
| 215 | IF (i /= i1) THEN |
---|
| 216 | IF ( (pttruetab(i1,Agrif_Procrank) /= pttruetab(i1,k)) .OR. & |
---|
| 217 | (cetruetab(i1,Agrif_Procrank) /= cetruetab(i1,k))) THEN |
---|
| 218 | tochangebis = .FALSE. |
---|
| 219 | EXIT |
---|
| 220 | ENDIF |
---|
| 221 | ENDIF |
---|
| 222 | ENDDO |
---|
[13027] | 223 | ! Strange CASE |
---|
| 224 | if ((pttruetab(i,k)>=pttruetab(i,Agrif_Procrank)).AND. & |
---|
| 225 | (cetruetab(i,k)<=cetruetab(i,Agrif_Procrank))) tochangebis = .FALSE. |
---|
| 226 | |
---|
| 227 | if (agrif_debug_interp) then |
---|
| 228 | print *,'tochangebis= ',tochangebis |
---|
| 229 | endif |
---|
[4777] | 230 | IF (tochangebis) THEN |
---|
| 231 | imin1 = max(pttruetab(i,Agrif_Procrank), pttruetab(i,k)) |
---|
| 232 | imax1 = min(cetruetab(i,Agrif_Procrank), cetruetab(i,k)) |
---|
| 233 | ! Always send the most interior points |
---|
[13027] | 234 | if (agrif_debug_interp) then |
---|
| 235 | print *,'imin1imax1= ',imin1,imax1 |
---|
| 236 | endif |
---|
[4777] | 237 | tochange = .false. |
---|
| 238 | IF (cetruetab(i,Agrif_Procrank) > cetruetab(i,k)) THEN |
---|
| 239 | DO j=imin1,imax1 |
---|
[13027] | 240 | IF ((bornesmax(i,k)-j) > (j-bornesmin(i,Agrif_Procrank))) THEN |
---|
[4777] | 241 | imintmp = j+1 |
---|
| 242 | tochange = .TRUE. |
---|
| 243 | ELSE |
---|
| 244 | EXIT |
---|
| 245 | ENDIF |
---|
| 246 | ENDDO |
---|
| 247 | ENDIF |
---|
| 248 | |
---|
| 249 | if (tochange) then |
---|
| 250 | pttruetab2(i,Agrif_Procrank) = imintmp |
---|
| 251 | endif |
---|
| 252 | |
---|
| 253 | tochange = .FALSE. |
---|
| 254 | imaxtmp=0 |
---|
| 255 | IF (pttruetab(i,Agrif_Procrank) < pttruetab(i,k)) THEN |
---|
| 256 | DO j=imax1,imin1,-1 |
---|
[13027] | 257 | IF ((j-bornesmin(i,k)) > (bornesmax(i,Agrif_Procrank)-j)) THEN |
---|
[4777] | 258 | imaxtmp = j-1 |
---|
| 259 | tochange = .TRUE. |
---|
| 260 | ELSE |
---|
| 261 | EXIT |
---|
| 262 | ENDIF |
---|
| 263 | ENDDO |
---|
| 264 | ENDIF |
---|
| 265 | |
---|
| 266 | if (tochange) then |
---|
| 267 | cetruetab2(i,Agrif_Procrank) = imaxtmp |
---|
| 268 | endif |
---|
| 269 | ENDIF |
---|
| 270 | enddo |
---|
| 271 | enddo |
---|
| 272 | |
---|
[13027] | 273 | if (agrif_debug_interp) then |
---|
| 274 | do k=0,Agrif_Nbprocs-1 |
---|
| 275 | print *,'Processeur ',k |
---|
| 276 | do i=1,nbdim |
---|
| 277 | print *,'ptcetretab2 = ',i,pttruetab2(i,k),cetruetab2(i,k) |
---|
| 278 | enddo |
---|
| 279 | enddo |
---|
| 280 | endif |
---|
| 281 | |
---|
[4777] | 282 | do k = 0,Agrif_NbProcs-1 |
---|
| 283 | ! |
---|
| 284 | sendtoproc(k) = .true. |
---|
| 285 | ! |
---|
[13027] | 286 | IF ( .not. memberoutall(k) ) THEN |
---|
| 287 | sendtoproc(k) = .false. |
---|
| 288 | ELSE |
---|
[4777] | 289 | !CDIR SHORTLOOP |
---|
| 290 | do i = 1,nbdim |
---|
| 291 | imin(i,k) = max(pttruetab2(i,Agrif_Procrank), pttruetabwhole(i,k)) |
---|
| 292 | imax(i,k) = min(cetruetab2(i,Agrif_Procrank), cetruetabwhole(i,k)) |
---|
| 293 | ! |
---|
| 294 | if ( (imin(i,k) > imax(i,k)) .and. (coords(i) /= 0) ) then |
---|
| 295 | sendtoproc(k) = .false. |
---|
| 296 | endif |
---|
| 297 | enddo |
---|
| 298 | ENDIF |
---|
| 299 | enddo |
---|
| 300 | ! |
---|
| 301 | call Exchangesamelevel_first(sendtoproc,nbdim,imin,imax,recvfromproc,imin_recv,imax_recv) |
---|
| 302 | !--------------------------------------------------------------------------------------------------- |
---|
| 303 | end subroutine Get_External_Data_first |
---|
| 304 | !=================================================================================================== |
---|
| 305 | ! |
---|
| 306 | !=================================================================================================== |
---|
| 307 | ! subroutine ExchangeSameLevel_first |
---|
| 308 | !--------------------------------------------------------------------------------------------------- |
---|
| 309 | subroutine ExchangeSameLevel_first ( sendtoproc, nbdim, imin, imax, recvfromproc, & |
---|
| 310 | imin_recv, imax_recv ) |
---|
| 311 | !--------------------------------------------------------------------------------------------------- |
---|
| 312 | LOGICAL, DIMENSION(0:Agrif_Nbprocs-1), intent(in) :: sendtoproc |
---|
| 313 | INTEGER, intent(in) :: nbdim |
---|
| 314 | INTEGER, DIMENSION(nbdim,0:Agrif_Nbprocs-1), intent(in) :: imin |
---|
| 315 | INTEGER, DIMENSION(nbdim,0:Agrif_Nbprocs-1), intent(in) :: imax |
---|
| 316 | LOGICAL, DIMENSION(0:Agrif_Nbprocs-1), intent(out) :: recvfromproc |
---|
| 317 | INTEGER, DIMENSION(nbdim,0:Agrif_Nbprocs-1), intent(out) :: imin_recv |
---|
| 318 | INTEGER, DIMENSION(nbdim,0:Agrif_Nbprocs-1), intent(out) :: imax_recv |
---|
| 319 | ! |
---|
| 320 | include 'mpif.h' |
---|
| 321 | INTEGER :: k |
---|
| 322 | INTEGER :: etiquette = 100 |
---|
| 323 | INTEGER :: code |
---|
| 324 | LOGICAL :: res |
---|
| 325 | INTEGER, DIMENSION(MPI_STATUS_SIZE) :: statut |
---|
| 326 | INTEGER, DIMENSION(nbdim,2,0:Agrif_Nbprocs-1) :: iminmax_temp |
---|
| 327 | |
---|
| 328 | do k = 0,Agrif_ProcRank-1 |
---|
| 329 | ! |
---|
| 330 | call MPI_SEND(sendtoproc(k),1,MPI_LOGICAL,k,etiquette,Agrif_mpi_comm,code) |
---|
| 331 | ! |
---|
| 332 | if (sendtoproc(k)) then |
---|
| 333 | iminmax_temp(:,1,k) = imin(:,k) |
---|
| 334 | iminmax_temp(:,2,k) = imax(:,k) |
---|
| 335 | call MPI_SEND(iminmax_temp(:,:,k),2*nbdim,MPI_INTEGER,k,etiquette,Agrif_mpi_comm,code) |
---|
| 336 | endif |
---|
| 337 | ! |
---|
| 338 | enddo |
---|
| 339 | ! |
---|
| 340 | ! Reception from others processors of the necessary part of the parent grid |
---|
| 341 | do k = Agrif_ProcRank+1,Agrif_Nbprocs-1 |
---|
| 342 | ! |
---|
| 343 | call MPI_RECV(res,1,MPI_LOGICAL,k,etiquette,Agrif_mpi_comm,statut,code) |
---|
| 344 | recvfromproc(k) = res |
---|
| 345 | ! |
---|
| 346 | if (recvfromproc(k)) then |
---|
| 347 | call MPI_RECV(iminmax_temp(:,:,k),2*nbdim,MPI_INTEGER,k,etiquette, & |
---|
| 348 | Agrif_mpi_comm,statut,code) |
---|
| 349 | imin_recv(:,k) = iminmax_temp(:,1,k) |
---|
| 350 | imax_recv(:,k) = iminmax_temp(:,2,k) |
---|
| 351 | endif |
---|
| 352 | ! |
---|
| 353 | enddo |
---|
| 354 | |
---|
| 355 | ! Reception from others processors of the necessary part of the parent grid |
---|
| 356 | do k = Agrif_ProcRank+1,Agrif_Nbprocs-1 |
---|
| 357 | ! |
---|
| 358 | call MPI_SEND(sendtoproc(k),1,MPI_LOGICAL,k,etiquette,Agrif_mpi_comm,code) |
---|
| 359 | ! |
---|
| 360 | if (sendtoproc(k)) then |
---|
| 361 | ! |
---|
| 362 | iminmax_temp(:,1,k) = imin(:,k) |
---|
| 363 | iminmax_temp(:,2,k) = imax(:,k) |
---|
| 364 | |
---|
| 365 | call MPI_SEND(iminmax_temp(:,:,k),2*nbdim,MPI_INTEGER,k,etiquette, & |
---|
| 366 | Agrif_mpi_comm,code) |
---|
| 367 | endif |
---|
| 368 | ! |
---|
| 369 | enddo |
---|
| 370 | ! |
---|
| 371 | ! |
---|
| 372 | ! Reception from others processors of the necessary part of the parent grid |
---|
| 373 | do k = Agrif_ProcRank-1,0,-1 |
---|
| 374 | ! |
---|
| 375 | call MPI_RECV(res,1,MPI_LOGICAL,k,etiquette,Agrif_mpi_comm,statut,code) |
---|
| 376 | recvfromproc(k) = res |
---|
| 377 | ! |
---|
| 378 | if (recvfromproc(k)) then |
---|
| 379 | ! |
---|
| 380 | call MPI_RECV(iminmax_temp(:,:,k),2*nbdim,MPI_INTEGER,k,etiquette, & |
---|
| 381 | Agrif_mpi_comm,statut,code) |
---|
| 382 | |
---|
| 383 | imin_recv(:,k) = iminmax_temp(:,1,k) |
---|
| 384 | imax_recv(:,k) = iminmax_temp(:,2,k) |
---|
| 385 | endif |
---|
| 386 | ! |
---|
| 387 | enddo |
---|
| 388 | !--------------------------------------------------------------------------------------------------- |
---|
| 389 | end subroutine ExchangeSamelevel_first |
---|
| 390 | !=================================================================================================== |
---|
| 391 | ! |
---|
| 392 | !=================================================================================================== |
---|
| 393 | ! subroutine ExchangeSameLevel |
---|
| 394 | !--------------------------------------------------------------------------------------------------- |
---|
| 395 | subroutine ExchangeSameLevel ( sendtoproc, recvfromproc, nbdim, & |
---|
| 396 | pttruetabwhole, cetruetabwhole, & |
---|
| 397 | imin, imax, imin_recv, imax_recv, & |
---|
| 398 | memberout, tempC, tempCextend ) |
---|
| 399 | !--------------------------------------------------------------------------------------------------- |
---|
| 400 | LOGICAL, DIMENSION(0:Agrif_Nbprocs-1), intent(in) :: sendtoproc |
---|
| 401 | LOGICAL, DIMENSION(0:Agrif_Nbprocs-1), intent(in) :: recvfromproc |
---|
| 402 | INTEGER, intent(in) :: nbdim |
---|
| 403 | INTEGER, DIMENSION(nbdim,0:Agrif_Nbprocs-1), intent(in) :: pttruetabwhole |
---|
| 404 | INTEGER, DIMENSION(nbdim,0:Agrif_Nbprocs-1), intent(in) :: cetruetabwhole |
---|
| 405 | INTEGER, DIMENSION(nbdim,0:Agrif_Nbprocs-1), intent(in) :: imin, imax |
---|
| 406 | INTEGER, DIMENSION(nbdim,0:Agrif_Nbprocs-1), intent(in) :: imin_recv, imax_recv |
---|
| 407 | LOGICAL, intent(in) :: memberout |
---|
| 408 | TYPE(Agrif_Variable), pointer, intent(inout) :: tempC, tempCextend |
---|
| 409 | ! |
---|
| 410 | include 'mpif.h' |
---|
| 411 | INTEGER :: i,k |
---|
| 412 | INTEGER :: etiquette = 100 |
---|
| 413 | INTEGER :: code, datasize |
---|
| 414 | INTEGER, DIMENSION(MPI_STATUS_SIZE) :: statut |
---|
| 415 | TYPE(Agrif_Variable), pointer, SAVE :: temprecv |
---|
| 416 | ! |
---|
| 417 | IF (memberout) THEN |
---|
| 418 | call Agrif_array_allocate(tempCextend, pttruetabwhole(:,Agrif_ProcRank), & |
---|
| 419 | cetruetabwhole(:,Agrif_ProcRank),nbdim) |
---|
| 420 | call Agrif_var_set_array_tozero(tempCextend,nbdim) |
---|
| 421 | ENDIF |
---|
| 422 | ! |
---|
[13027] | 423 | if (agrif_debug_interp) then |
---|
| 424 | print *,'PROCESSEUR = ',Agrif_Procrank |
---|
| 425 | print *,'SENDTOPROC = ',sendtoproc(Agrif_Procrank) |
---|
| 426 | if (sendtoproc(Agrif_Procrank)) then |
---|
| 427 | print *,'imin imax = ',imin(:,Agrif_Procrank),imax(:,Agrif_Procrank) |
---|
| 428 | endif |
---|
| 429 | endif |
---|
[4777] | 430 | IF (sendtoproc(Agrif_ProcRank)) THEN |
---|
| 431 | call Agrif_var_copy_array(tempCextend,imin(:,Agrif_Procrank),imax(:,Agrif_Procrank), & |
---|
| 432 | tempC, imin(:,Agrif_Procrank),imax(:,Agrif_Procrank), & |
---|
| 433 | nbdim) |
---|
| 434 | ENDIF |
---|
| 435 | ! |
---|
| 436 | do k = 0,Agrif_ProcRank-1 |
---|
| 437 | ! |
---|
| 438 | if (sendtoproc(k)) then |
---|
| 439 | ! |
---|
| 440 | datasize = 1 |
---|
| 441 | ! |
---|
| 442 | !CDIR SHORTLOOP |
---|
| 443 | do i = 1,nbdim |
---|
| 444 | datasize = datasize * (imax(i,k)-imin(i,k)+1) |
---|
| 445 | enddo |
---|
| 446 | ! |
---|
| 447 | SELECT CASE(nbdim) |
---|
| 448 | CASE(1) |
---|
| 449 | call MPI_SEND(tempC%array1(imin(1,k):imax(1,k)), & |
---|
| 450 | datasize,Agrif_MPI_prec,k,etiquette, & |
---|
| 451 | Agrif_mpi_comm,code) |
---|
| 452 | CASE(2) |
---|
| 453 | call MPI_SEND(tempC%array2(imin(1,k):imax(1,k), & |
---|
| 454 | imin(2,k):imax(2,k)), & |
---|
| 455 | datasize,Agrif_MPI_prec,k,etiquette, & |
---|
| 456 | Agrif_mpi_comm,code) |
---|
| 457 | CASE(3) |
---|
| 458 | call Agrif_Send_3Darray(tempC%array3,lbound(tempC%array3),imin(:,k),imax(:,k),k) |
---|
| 459 | CASE(4) |
---|
| 460 | call MPI_SEND(tempC%array4(imin(1,k):imax(1,k), & |
---|
| 461 | imin(2,k):imax(2,k), & |
---|
| 462 | imin(3,k):imax(3,k), & |
---|
| 463 | imin(4,k):imax(4,k)), & |
---|
| 464 | datasize,Agrif_MPI_prec,k,etiquette, & |
---|
| 465 | Agrif_mpi_comm,code) |
---|
| 466 | CASE(5) |
---|
| 467 | call MPI_SEND(tempC%array5(imin(1,k):imax(1,k), & |
---|
| 468 | imin(2,k):imax(2,k), & |
---|
| 469 | imin(3,k):imax(3,k), & |
---|
| 470 | imin(4,k):imax(4,k), & |
---|
| 471 | imin(5,k):imax(5,k)), & |
---|
| 472 | datasize,Agrif_MPI_prec,k,etiquette, & |
---|
| 473 | Agrif_mpi_comm,code) |
---|
| 474 | CASE(6) |
---|
| 475 | call MPI_SEND(tempC%array6(imin(1,k):imax(1,k), & |
---|
| 476 | imin(2,k):imax(2,k), & |
---|
| 477 | imin(3,k):imax(3,k), & |
---|
| 478 | imin(4,k):imax(4,k), & |
---|
| 479 | imin(5,k):imax(5,k), & |
---|
| 480 | imin(6,k):imax(6,k)), & |
---|
| 481 | datasize,Agrif_MPI_prec,k,etiquette, & |
---|
| 482 | Agrif_mpi_comm,code) |
---|
| 483 | END SELECT |
---|
| 484 | ! |
---|
| 485 | endif |
---|
| 486 | enddo |
---|
| 487 | ! |
---|
| 488 | ! Reception from others processors of the necessary part of the parent grid |
---|
| 489 | do k = Agrif_ProcRank+1,Agrif_Nbprocs-1 |
---|
| 490 | ! |
---|
| 491 | if (recvfromproc(k)) then |
---|
| 492 | ! |
---|
| 493 | datasize = 1 |
---|
| 494 | ! |
---|
| 495 | !CDIR SHORTLOOP |
---|
| 496 | do i = 1,nbdim |
---|
| 497 | datasize = datasize * (imax_recv(i,k)-imin_recv(i,k)+1) |
---|
| 498 | enddo |
---|
| 499 | |
---|
| 500 | if (.not.associated(temprecv)) allocate(temprecv) |
---|
| 501 | call Agrif_array_allocate(temprecv,imin_recv(:,k),imax_recv(:,k),nbdim) |
---|
| 502 | |
---|
| 503 | SELECT CASE(nbdim) |
---|
| 504 | CASE(1) |
---|
| 505 | call MPI_RECV(temprecv%array1,datasize,Agrif_MPI_prec,k,etiquette, & |
---|
| 506 | Agrif_mpi_comm,statut,code) |
---|
| 507 | CASE(2) |
---|
| 508 | call MPI_RECV(temprecv%array2,datasize,Agrif_MPI_prec,k,etiquette, & |
---|
| 509 | Agrif_mpi_comm,statut,code) |
---|
| 510 | CASE(3) |
---|
| 511 | call MPI_RECV(temprecv%array3,datasize,Agrif_MPI_prec,k,etiquette, & |
---|
| 512 | Agrif_mpi_comm,statut,code) |
---|
| 513 | CASE(4) |
---|
| 514 | call MPI_RECV(temprecv%array4,datasize,Agrif_MPI_prec,k,etiquette, & |
---|
| 515 | Agrif_mpi_comm,statut,code) |
---|
| 516 | CASE(5) |
---|
| 517 | call MPI_RECV(temprecv%array5,datasize,Agrif_MPI_prec,k,etiquette, & |
---|
| 518 | Agrif_mpi_comm,statut,code) |
---|
| 519 | CASE(6) |
---|
| 520 | call MPI_RECV(temprecv%array6,datasize,Agrif_MPI_prec,k,etiquette, & |
---|
| 521 | Agrif_mpi_comm,statut,code) |
---|
| 522 | END SELECT |
---|
| 523 | |
---|
| 524 | call Agrif_var_replace_value(tempCextend,temprecv,imin_recv(:,k),imax_recv(:,k),0.,nbdim) |
---|
| 525 | call Agrif_array_deallocate(temprecv,nbdim) |
---|
| 526 | ! |
---|
| 527 | endif |
---|
| 528 | enddo |
---|
| 529 | |
---|
| 530 | ! Reception from others processors of the necessary part of the parent grid |
---|
| 531 | do k = Agrif_ProcRank+1,Agrif_Nbprocs-1 |
---|
| 532 | ! |
---|
| 533 | if (sendtoproc(k)) then |
---|
| 534 | ! |
---|
| 535 | SELECT CASE(nbdim) |
---|
| 536 | CASE(1) |
---|
| 537 | datasize=SIZE(tempC%array1(imin(1,k):imax(1,k))) |
---|
| 538 | call MPI_SEND(tempC%array1(imin(1,k):imax(1,k)), & |
---|
| 539 | datasize,Agrif_MPI_prec,k,etiquette, & |
---|
| 540 | Agrif_mpi_comm,code) |
---|
| 541 | CASE(2) |
---|
| 542 | datasize=SIZE(tempC%array2(imin(1,k):imax(1,k), & |
---|
| 543 | imin(2,k):imax(2,k))) |
---|
| 544 | call MPI_SEND(tempC%array2(imin(1,k):imax(1,k), & |
---|
| 545 | imin(2,k):imax(2,k)), & |
---|
| 546 | datasize,Agrif_MPI_prec,k,etiquette, & |
---|
| 547 | Agrif_mpi_comm,code) |
---|
| 548 | CASE(3) |
---|
| 549 | datasize=SIZE(tempC%array3(imin(1,k):imax(1,k), & |
---|
| 550 | imin(2,k):imax(2,k), & |
---|
| 551 | imin(3,k):imax(3,k))) |
---|
| 552 | call MPI_SEND(tempC%array3(imin(1,k):imax(1,k), & |
---|
| 553 | imin(2,k):imax(2,k), & |
---|
| 554 | imin(3,k):imax(3,k)), & |
---|
| 555 | datasize,Agrif_MPI_prec,k,etiquette, & |
---|
| 556 | Agrif_mpi_comm,code) |
---|
| 557 | CASE(4) |
---|
| 558 | datasize=SIZE(tempC%array4(imin(1,k):imax(1,k), & |
---|
| 559 | imin(2,k):imax(2,k), & |
---|
| 560 | imin(3,k):imax(3,k), & |
---|
| 561 | imin(4,k):imax(4,k))) |
---|
| 562 | call MPI_SEND(tempC%array4(imin(1,k):imax(1,k), & |
---|
| 563 | imin(2,k):imax(2,k), & |
---|
| 564 | imin(3,k):imax(3,k), & |
---|
| 565 | imin(4,k):imax(4,k)), & |
---|
| 566 | datasize,Agrif_MPI_prec,k,etiquette, & |
---|
| 567 | Agrif_mpi_comm,code) |
---|
| 568 | CASE(5) |
---|
| 569 | datasize=SIZE(tempC%array5(imin(1,k):imax(1,k), & |
---|
| 570 | imin(2,k):imax(2,k), & |
---|
| 571 | imin(3,k):imax(3,k), & |
---|
| 572 | imin(4,k):imax(4,k), & |
---|
| 573 | imin(5,k):imax(5,k))) |
---|
| 574 | call MPI_SEND(tempC%array5(imin(1,k):imax(1,k), & |
---|
| 575 | imin(2,k):imax(2,k), & |
---|
| 576 | imin(3,k):imax(3,k), & |
---|
| 577 | imin(4,k):imax(4,k), & |
---|
| 578 | imin(5,k):imax(5,k)), & |
---|
| 579 | datasize,Agrif_MPI_prec,k,etiquette, & |
---|
| 580 | Agrif_mpi_comm,code) |
---|
| 581 | CASE(6) |
---|
| 582 | datasize=SIZE(tempC%array6(imin(1,k):imax(1,k), & |
---|
| 583 | imin(2,k):imax(2,k), & |
---|
| 584 | imin(3,k):imax(3,k), & |
---|
| 585 | imin(4,k):imax(4,k), & |
---|
| 586 | imin(5,k):imax(5,k), & |
---|
| 587 | imin(6,k):imax(6,k))) |
---|
| 588 | call MPI_SEND(tempC%array6(imin(1,k):imax(1,k), & |
---|
| 589 | imin(2,k):imax(2,k), & |
---|
| 590 | imin(3,k):imax(3,k), & |
---|
| 591 | imin(4,k):imax(4,k), & |
---|
| 592 | imin(5,k):imax(5,k), & |
---|
| 593 | imin(6,k):imax(6,k)), & |
---|
| 594 | datasize,Agrif_MPI_prec,k,etiquette, & |
---|
| 595 | Agrif_mpi_comm,code) |
---|
| 596 | END SELECT |
---|
| 597 | ! |
---|
| 598 | endif |
---|
| 599 | ! |
---|
| 600 | enddo |
---|
| 601 | ! |
---|
| 602 | ! Reception from others processors of the necessary part of the parent grid |
---|
| 603 | do k = Agrif_ProcRank-1,0,-1 |
---|
| 604 | ! |
---|
| 605 | if (recvfromproc(k)) then |
---|
| 606 | ! |
---|
| 607 | if (.not.associated(temprecv)) allocate(temprecv) |
---|
| 608 | call Agrif_array_allocate(temprecv,imin_recv(:,k),imax_recv(:,k),nbdim) |
---|
| 609 | |
---|
| 610 | SELECT CASE(nbdim) |
---|
| 611 | CASE(1) |
---|
| 612 | datasize=SIZE(temprecv%array1) |
---|
| 613 | call MPI_RECV(temprecv%array1,datasize,Agrif_MPI_prec,k,etiquette,& |
---|
| 614 | Agrif_mpi_comm,statut,code) |
---|
| 615 | CASE(2) |
---|
| 616 | datasize=SIZE(temprecv%array2) |
---|
| 617 | call MPI_RECV(temprecv%array2,datasize,Agrif_MPI_prec,k,etiquette,& |
---|
| 618 | Agrif_mpi_comm,statut,code) |
---|
| 619 | CASE(3) |
---|
| 620 | datasize=SIZE(temprecv%array3) |
---|
| 621 | call MPI_RECV(temprecv%array3,datasize,Agrif_MPI_prec,k,etiquette,& |
---|
| 622 | Agrif_mpi_comm,statut,code) |
---|
| 623 | CASE(4) |
---|
| 624 | datasize=SIZE(temprecv%array4) |
---|
| 625 | call MPI_RECV(temprecv%array4,datasize,Agrif_MPI_prec,k,etiquette,& |
---|
| 626 | Agrif_mpi_comm,statut,code) |
---|
| 627 | CASE(5) |
---|
| 628 | datasize=SIZE(temprecv%array5) |
---|
| 629 | call MPI_RECV(temprecv%array5,datasize,Agrif_MPI_prec,k,etiquette,& |
---|
| 630 | Agrif_mpi_comm,statut,code) |
---|
| 631 | CASE(6) |
---|
| 632 | datasize=SIZE(temprecv%array6) |
---|
| 633 | call MPI_RECV(temprecv%array6,datasize,Agrif_MPI_prec,k,etiquette,& |
---|
| 634 | Agrif_mpi_comm,statut,code) |
---|
| 635 | END SELECT |
---|
| 636 | |
---|
| 637 | call Agrif_var_replace_value(tempCextend,temprecv,imin_recv(:,k),imax_recv(:,k),0.,nbdim) |
---|
| 638 | call Agrif_array_deallocate(temprecv,nbdim) |
---|
| 639 | ! |
---|
| 640 | endif |
---|
| 641 | ! |
---|
| 642 | enddo |
---|
| 643 | !--------------------------------------------------------------------------------------------------- |
---|
| 644 | end subroutine ExchangeSamelevel |
---|
| 645 | !=================================================================================================== |
---|
| 646 | ! |
---|
| 647 | !=================================================================================================== |
---|
| 648 | ! subroutine Agrif_Send_3Darray |
---|
| 649 | !--------------------------------------------------------------------------------------------------- |
---|
| 650 | subroutine Agrif_Send_3Darray ( tab3D, bounds, imin, imax, k ) |
---|
| 651 | !--------------------------------------------------------------------------------------------------- |
---|
| 652 | integer, dimension(3), intent(in) :: bounds |
---|
| 653 | real, dimension(bounds(1):,bounds(2):,bounds(3):), target, intent(in) :: tab3D |
---|
| 654 | integer, dimension(3), intent(in) :: imin, imax |
---|
| 655 | integer, intent(in) :: k |
---|
| 656 | ! |
---|
| 657 | integer :: etiquette = 100 |
---|
| 658 | integer :: datasize, code |
---|
| 659 | include 'mpif.h' |
---|
| 660 | |
---|
| 661 | datasize = SIZE(tab3D(imin(1):imax(1), & |
---|
| 662 | imin(2):imax(2), & |
---|
| 663 | imin(3):imax(3))) |
---|
| 664 | |
---|
| 665 | call MPI_SEND( tab3D( imin(1):imax(1), & |
---|
| 666 | imin(2):imax(2), & |
---|
| 667 | imin(3):imax(3)), & |
---|
| 668 | datasize,Agrif_MPI_prec,k,etiquette,Agrif_mpi_comm,code) |
---|
| 669 | !--------------------------------------------------------------------------------------------------- |
---|
| 670 | end subroutine Agrif_Send_3Darray |
---|
| 671 | !=================================================================================================== |
---|
| 672 | ! |
---|
| 673 | #else |
---|
| 674 | subroutine dummy_Agrif_Mpp () |
---|
| 675 | end subroutine dummy_Agrif_Mpp |
---|
| 676 | #endif |
---|
| 677 | ! |
---|
| 678 | end Module Agrif_Mpp |
---|